include file fortran

Barry Smith bsmith at mcs.anl.gov
Mon May 18 08:23:59 CDT 2009


    The calling sequence was changed for VecSetOption and MatSetOption  
see http://www.mcs.anl.gov/petsc/petsc-as/documentation/changes/300.html

    Barry

On May 18, 2009, at 7:00 AM, Tahar Amari wrote:

> Hello,
>
> I finally succeed building Petsc and the fortran  code thanks to  
> your advice .
>
> However when I run the  code (piece of it is below)
>
>     	if(pe==0) write(*,*) "before MatCopy"
>      CALL MatCopy(dvoeta_petsc,lhs_petsc,SAME_NONZERO_PATTERN,
>     & ierror)
>     	if(pe==0) write(*,*) "after MatCopy"
>
>
> I got the following message (in between my printing for debugging) .
>
> before MatCopy
>
> [0]PETSC ERROR: --------------------- Error Message  
> ------------------------------------
> [0]PETSC ERROR: Corrupt argument: see http://www.mcs.anl.gov/petsc/petsc-as/documentation/troubleshooting.html#Corrupt 
> !
> [0]PETSC ERROR: Invalid Pointer to Object: Parameter # 1!
> [0]PETSC ERROR:  
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.0.0, Patch 5, Mon Apr 13  
> 09:15:37 CDT 2009
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR:  
> ------------------------------------------------------------------------
> [0]PETSC ERROR: mh4d on a macx named imac-de-tahar-amari.local by  
> amari Mon May 18 13:47:42 2009
> [0]PETSC ERROR: Libraries linked from /Data/Poub1/petsc-3.0.0-p5/ 
> macx/lib
> [0]PETSC ERROR: Configure run at Sat May 16 13:49:21 2009
> [0]PETSC ERROR: Configure options -with-cc=gcc --download-mpich=1 -- 
> download-f-bl
>
> after MatCopy
>
>
> I looked at the piece of code which define the first argument  
> "dvoeta_petsc"
>
> which is
>
>         Mat, SAVE :: dvoeta_petsc
>          CALL init_petsc_mat_vector(dvoeta_petsc,diagonal=dvoeta)
>
>
> This code was working with Petsc version 2.xxx
> I wonder if some changes to version 3 might  be the reason ?
>
>
> May it be the call to
>
>        CALL MatSetOption(matrix,MAT_NO_NEW_NONZERO_LOCATIONS,ierr)
>
> in this function please ?
>
> c     subprogram 8. init_petsc_mat_vector
> c     Make a petsc matrix out of an operator
> c 
> -----------------------------------------------------------------------
>      SUBROUTINE init_petsc_mat_vector(matrix,v2v_operator_vv,diagonal,
>     &bctype)
> c 
> -----------------------------------------------------------------------
>      USE operator_mod
> #include "include/finclude/petsc.h"
> #include "include/finclude/petscvec.h"
> #include "include/finclude/petscmat.h"
> #include "include/finclude/petscao.h"
>
> c 
> -----------------------------------------------------------------------
>      Mat, INTENT(INOUT) :: matrix
>      INTERFACE
>        SUBROUTINE v2v_operator_vv(x)
>          USE operator_mod
>          TYPE (v2v_operators), INTENT(INOUT) :: x(:)
>        END SUBROUTINE v2v_operator_vv
>      END INTERFACE
>      OPTIONAL :: v2v_operator_vv
>      TYPE (vertex_scalar), INTENT(IN), OPTIONAL :: diagonal
>      INTEGER(i4), INTENT(IN), OPTIONAL :: bctype
>      INTEGER(i4) :: m,n,ierr,i,n1,m1,nv_local_inside,
>     &  iv_local,iv_couple,d_nz,o_nz
>      TYPE(v2v_operators),  ALLOCATABLE :: v2voperator(:)
>      REAL(r8) :: v(3,3)
>      LOGICAL, SAVE :: first=.true.
> c 
> -----------------------------------------------------------------------
>
>      CALL get_grid(thegrid,nv_local_inside=nv_local_inside)
>      CALL get_nonzero(thegrid,d_nz,o_nz)
>
>      CALL MatCreateMPIAIJ(comm_grid,3*nv_local_inside,
>     & 3*nv_local_inside,PETSC_DECIDE,PETSC_DECIDE,3*d_nz,
>     & PETSC_NULL_INTEGER,3*o_nz,PETSC_NULL_INTEGER,matrix,ierr)
>
>      ALLOCATE(v2voperator(nv_local_inside))
> c 
> -----------------------------------------------------------------------
> c See if it is a diagonal operator or not.
> c 
> -----------------------------------------------------------------------
>      IF (PRESENT(diagonal)) THEN
>        IF(first) CALL terminator("Error in init_petsc_mat_vector: "//
>     &  "You shouldn't call it the first time with a diagonal  
> operator")
>        IF(PRESENT(bctype)) THEN
>          CALL v2v_diagonal_vv(v2voperator,diagonal,bctype)
>        ELSE
>          CALL v2v_diagonal_vv(v2voperator,diagonal,1)
>        ENDIF
>      ELSE
>        CALL v2v_operator_vv(v2voperator)
>      ENDIF
>
>      DO iv_local=1,nv_local_inside
>        CALL index_vv2petsc(m,iv_local,1)
>        DO i=1,SIZE(v2voperator(iv_local)%neighbors)
>          iv_couple=v2voperator(iv_local)%neighbors(i)
>          CALL index_vv2petsc(n,iv_couple,1)
>          v=v2voperator(iv_local)%coupling_with(i)%a
>          DO m1=0,2
>            DO n1=0,2
>              CALL MatSetValues(matrix,1,m+m1,1,n+n1,v(1+m1,1+n1),
>     &        INSERT_VALUES,ierr)
>            ENDDO
>          ENDDO
>        ENDDO
>      ENDDO
>
>      CALL MatAssemblyBegin(matrix,MAT_FINAL_ASSEMBLY,ierr)
>      CALL MatAssemblyEnd(matrix,MAT_FINAL_ASSEMBLY,ierr)
>
>      IF(first) THEN
>        CALL MatSetOption(matrix,MAT_NO_NEW_NONZERO_LOCATIONS,ierr)
>        first=.false.
>      ENDIF
>
>      CALL destroy(v2voperator)
>      DEALLOCATE(v2voperator)
>
>      RETURN
>
>
>
>



More information about the petsc-users mailing list