[petsc-users] Solver compilation with 64-bit version of PETSc under Windows 10 using Cygwin

Smith, Barry F. bsmith at mcs.anl.gov
Fri Jan 17 15:55:05 CST 2020



   I am working on it. It requires a large number of fixes; I hope to have it working by tonight.

   Barry


> On Jan 17, 2020, at 2:40 AM, Дмитрий Мельничук <dmitry.melnichuk at geosteertech.com> wrote:
> 
> Thank you for your replies.
> 
> Tried to configure commited version of PETSc from Satish Balay branch (balay/fix-ftn-i8/maint) and faced to the same error when running test example ex5f
>  
> 
> call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
>  
> Error: Type mismatch in argument «z» at (1); passed INTEGER(4) to INTEGER(8)
>  
> 
> At the moment, some subroutines (such as PetscInitialize, PetscFinalize, MatSetValue, VecSetValue) work with the correct size of the variable ierr defined as PetscErrorCode, and some do not.
> The following subroutines still require ierr to be of type INTEGER(8):
>  
> 
> VecGetSubVector, VecAssemblyBegin, VecAssemblyEnd, VecScatterBegin, VecScatterEnd, VecScatterDestroy, VecCreateMPI, VecDuplicate, VecZeroEntries, VecAYPX, VecWAXPY, VecWAXPY
> MatMult, MatDestroy, MatAssemblyBegin, MatAssemblyEnd, MatZeroEntries, MatCreateSubMatrix, MatScale, MatDiagonalSet, MatGetDiagonal, MatDuplicat, MatSetSizes, MatSetFromOptions
> 
> Unfortunately, I'm not sure if this is the only issue that occurs when switching to the 64-bit version of PETSc.
> I can set the size of the variables ierr so that the solver compilation process completes successfully, but I get the following error when solving linear algebra system by use of KSPSolve subroutine:
>  
>  
> [0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range
>  
> Because the solver with the 32-bit version of PETSc works properly, I suppose that the cause of the errors (for 64-bit version of PETSc) is the inappropriate size of the variables.
> I compiled PETSc with flags: with-64-bit-indices and -fdefault-integer-8. 
> Also changed the size of MPI_Integer to MPI_Integer8:
> MPI_Bcast(npart,nnds,MPI_Integer8,0,MPI_Comm_World,ierr).
>  
> I am probably missing something else.
> 
> ​​​​​​​
> Kind regards,
> Dmitry Melnichuk
>  
> 16.01.2020, 01:26, "Balay, Satish" <balay at mcs.anl.gov>:
> I have some changes (incomplete) here -
> 
> my hack to bfort.
> 
> diff --git a/src/bfort/bfort.c b/src/bfort/bfort.c
> index 0efe900..31ff154 100644
> --- a/src/bfort/bfort.c
> +++ b/src/bfort/bfort.c
> @@ -1654,7 +1654,7 @@ void PrintDefinition( FILE *fout, int is_function, char *name, int nstrings,
> 
>      /* Add a "decl/result(name) for functions */
>      if (useFerr) {
> - OutputFortranToken( fout, 7, "integer" );
> + OutputFortranToken( fout, 7, "PetscErrorCode" );
>         OutputFortranToken( fout, 1, errArgNameParm);
>      } else if (is_function) {
>         OutputFortranToken( fout, 7, ArgToFortran( rt->name ) );
> 
> 
> And my changes to petsc are on branch balay/fix-ftn-i8/maint
> 
> Satish
> 
> On Wed, 15 Jan 2020, Smith, Barry F. via petsc-users wrote:
>  
> 
> 
>    Working on it now; may be doable
> 
> 
> 
>  > On Jan 15, 2020, at 11:55 AM, Matthew Knepley <knepley at gmail.com> wrote:
>  >
>  > On Wed, Jan 15, 2020 at 10:26 AM Дмитрий Мельничук <dmitry.melnichuk at geosteertech.com> wrote:
>  > > And I'm not sure why you are having to use PetscInt for ierr. All PETSc routines should be suing 'PetscErrorCode for ierr'
>  >
>  > If I define ierr as PetscErrorCode for all subroutines given below
>  >
>  > call VecDuplicate(Vec_U,Vec_Um,ierr)
>  > call VecCopy(Vec_U,Vec_Um,ierr)
>  > call VecGetLocalSize(Vec_U,j,ierr)
>  > call VecGetOwnershipRange(Vec_U,j1,j2,ierr)
>  >
>  > then errors occur with first three subroutines:
>  > Error: Type mismatch in argument «z» at (1); passed INTEGER(4) to INTEGER(8).
>  >
>  > Barry,
>  >
>  > It looks like the ftn-auto interfaces are using 'integer' for the error code, whereas the ftn-custom is using PetscErrorCode.
>  > Could we make the generated ones use integer?
>  >
>  > Thanks,
>  >
>  > Matt
>  >
>  > Therefore I was forced to define ierr as PetscInt for VecDuplicate, VecCopy, VecGetLocalSize subroutines to fix these errors.
>  > Why some subroutines sue 8-bytes integer type of ierr (PetscInt), while others - 4-bytes integer type of ierr (PetscErrorCode) remains a mystery for me.
>  >
>  > > What version of PETSc are you using?
>  >
>  > version 3.12.2
>  >
>  > > Are you seeing this issue with a PETSc example?
>  >
>  > I will check it tomorrow and let you know.
>  >
>  > Kind regards,
>  > Dmitry Melnichuk
>  >
>  >
>  >
>  > 15.01.2020, 17:14, "Balay, Satish" <balay at mcs.anl.gov>:
>  > -fdefault-integer-8 is likely to break things [esp with MPI - where 'integer' is used everywhere for ex - MPI_Comm etc - so MPI includes become incompatible with the MPI library with -fdefault-integer-8.]
>  >
>  > And I'm not sure why you are having to use PetscInt for ierr. All PETSc routines should be suing 'PetscErrorCode for ierr'
>  >
>  > What version of PETSc are you using? Are you seeing this issue with a PETSc example?
>  >
>  > Satish
>  >
>  > On Wed, 15 Jan 2020, Дмитрий Мельничук wrote:
>  >
>  >
>  > Hello all!
>  > At present time I need to compile solver called Defmod (https://bitbucket.org/stali/defmod/wiki/Home), which is written in Fortran 95.
>  > Defmod uses PETSc for solving linear algebra system.
>  > Solver compilation with 32-bit version of PETSc does not cause any problem.
>  > But solver compilation with 64-bit version of PETSc produces an error with size of ierr PETSc variable.
>  >
>  > 1. For example, consider the following statements written in Fortran:
>  >
>  >
>  > PetscErrorCode :: ierr_m
>  > PetscInt :: ierr
>  > ...
>  > ...
>  > call VecDuplicate(Vec_U,Vec_Um,ierr)
>  > call VecCopy(Vec_U,Vec_Um,ierr)
>  > call VecGetLocalSize(Vec_U,j,ierr)
>  > call VecGetOwnershipRange(Vec_U,j1,j2,ierr_m)
>  >
>  >
>  > As can be seen first three subroutunes require ierr to be size of INTEGER(8), while the last subroutine (VecGetOwnershipRange) requires ierr to be size of INTEGER(4).
>  > Using the same integer format gives an error:
>  >
>  > There is no specific subroutine for the generic ‘vecgetownershiprange’ at (1)
>  >
>  > 2. Another example is:
>  >
>  >
>  > call MatAssemblyBegin(Mat_K,Mat_Final_Assembly,ierr)
>  > CHKERRA(ierr)
>  > call MatAssemblyEnd(Mat_K,Mat_Final_Assembly,ierr)
>  >
>  >
>  > I am not able to define an appropriate size if ierr in CHKERRA(ierr). If I choose INTEGER(8), the error "Type mismatch in argument ‘ierr’ at (1); passed INTEGER(8) to
>  > INTEGER(4)" occurs.
>  > If I define ierr as INTEGER(4), the error "Type mismatch in argument ‘ierr’ at (1); passed INTEGER(4) to INTEGER(8)" appears.
>  >
>  > 3. If I change the sizes of ierr vaiables as error messages require, the compilation completed successfully, but an error occurs when calculating the RHS vector with
>  > following message:
>  > [0]PETSC ERROR: Out of range index value -4 cannot be negative
>  >
>  >
>  > Command to configure 32-bit version of PETSc under Windows 10 using Cygwin:
>  > ./configure --with-cc=x86_64-w64-mingw32-gcc --with-cxx=x86_64-w64-mingw32-g++ --with-fc=x86_64-w64-mingw32-gfortran --download-fblaslapack
>  > --with-mpi-include=/cygdrive/c/MPISDK/Include --with-mpi-lib=/cygdrive/c/MPISDK/Lib/libmsmpi.a --with-mpi-mpiexec=/cygdrive/c/MPI/Bin/mpiexec.exe --with-debugging=yes
>  > -CFLAGS='-O2' -CXXFLAGS='-O2' -FFLAGS='-O2 -static-libgfortran -static -lpthread -fno-range-check' --with-shared-libraries=no
>  > Command to configure 64-bit version of PETSc under Windows 10 using Cygwin:./configure --with-cc=x86_64-w64-mingw32-gcc --with-cxx=x86_64-w64-mingw32-g++
>  > --with-fc=x86_64-w64-mingw32-gfortran --download-fblaslapack --with-mpi-include=/cygdrive/c/MPISDK/Include --with-mpi-lib=/cygdrive/c/MPISDK/Lib/libmsmpi.a
>  > --with-mpi-mpiexec=/cygdrive/c/MPI/Bin/mpiexec.exe --with-debugging=yes -CFLAGS='-O2' -CXXFLAGS='-O2' -FFLAGS='-O2 -static-libgfortran -static -lpthread -fno-range-check
>  > -fdefault-integer-8' --with-shared-libraries=no --with-64-bit-indices --known-64-bit-blas-indices
>  >
>  >
>  > Kind regards,
>  > Dmitry Melnichuk
>  >
>  >
>  >
>  >
>  > --
>  > What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>  > -- Norbert Wiener
>  >
>  > https://www.cse.buffalo.edu/~knepley/
> 
>  



More information about the petsc-users mailing list