[petsc-users] Problem with SNESsetFunction()

Barry Smith bsmith at petsc.dev
Sat Oct 31 09:14:02 CDT 2020



> On Oct 30, 2020, at 10:11 PM, baikadi pranay <pranayreddy865 at gmail.com> wrote:
> 
> Hello,
> I have a couple of questions regarding SNESSetFunction usage, when programming in Fortran90.
> 
> 1) I have the following usage paradigm.
>    call SNESSetFunction(snes,f_non,FormFunction,0,ierr)
>    subroutine FormFunction(snes,x,r,dummy,ierr)
> In the FormFunction subroutine, the function values are stored in the vector r. I see that these values are formed correctly. But when I use FormFunction in  SNESSetFunction(), the values are not getting populated into f_non and all of the values in f_non are zero. 

> Should the name of the variable used to store the function value be same in  SNESSetFunction and FormFunction?

   It does not need to be the same, they are just the variables in each function

> And should I be calling the SNESComputeFunction() after calling SNESSetFunction()?

   No, that is a developer function called in PETSc, one would not need to call that.

> 
> 2) In the subroutine FormFunction, should the vector objects created be destroyed before ending the subroutine?

   What vectors? If you are creating any work vectors you need within the the FormFunction, yes those should be destroyed. But not the input and output functions.

   Here is any example from src/snes/tutorials/ex1f.F90 Note you call VecGetArrayF90() to access the arrays for the vectors, put the values into the arrays


subroutine FormFunction(snes,X,F,user,ierr)
      implicit none

!  Input/output variables:
      SNES           snes
      Vec            X,F
      PetscErrorCode ierr
      type (userctx) user
      DM             da

!  Declarations for use with local arrays:
      PetscScalar,pointer :: lx_v(:),lf_v(:)
      Vec            localX

!  Scatter ghost points to local vector, using the 2-step process
!     DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
!  By placing code between these two statements, computations can
!  be done while messages are in transition.
      call SNESGetDM(snes,da,ierr);CHKERRQ(ierr)
      call DMGetLocalVector(da,localX,ierr);CHKERRQ(ierr)
      call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr);CHKERRQ(ierr)
      call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr);CHKERRQ(ierr)

!  Get a pointer to vector data.
!    - For default PETSc vectors, VecGetArray90() returns a pointer to
!      the data array. Otherwise, the routine is implementation dependent.
!    - You MUST call VecRestoreArrayF90() when you no longer need access to
!      the array.
!    - Note that the interface to VecGetArrayF90() differs from VecGetArray(),
!      and is useable from Fortran-90 Only.

      call VecGetArrayF90(localX,lx_v,ierr);CHKERRQ(ierr)
      call VecGetArrayF90(F,lf_v,ierr);CHKERRQ(ierr)

!  Compute function over the locally owned part of the grid
      call FormFunctionLocal(lx_v,lf_v,user,ierr);CHKERRQ(ierr)

!  Restore vectors
      call VecRestoreArrayF90(localX,lx_v,ierr);CHKERRQ(ierr)
      call VecRestoreArrayF90(F,lf_v,ierr);CHKERRQ(ierr)

!  Insert values into global vector

      call DMRestoreLocalVector(da,localX,ierr);CHKERRQ(ierr)
      call PetscLogFlops(11.0d0*user%ym*user%xm,ierr)

!      call VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)
!      call VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)
      return
      end subroutine formfunction
      end module f90module



> 
> Please let me know if you need any further information. Thank you in advance.
> Best regards,
> Pranay.
> 
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20201031/cea604b3/attachment.html>


More information about the petsc-users mailing list