Jose,
   I've tried lots of combinations but I still get the same error. I think the signatures are all correct.  I've attached the routine in case you see something obvious.  If not I will try to make a standalone program that generates the same compile error.

   upremas.F:66:72:

       66 |      &                               ierr)
   | 1
   Error: There is no specific subroutine for the generic
   'matmpibaijsetpreallocation' at (1)
   upremas.F:69:72:

       69 |      &                               ierr)
   | 1
   Error: There is no specific subroutine for the generic
   'matseqbaijsetpreallocation' at (1)
   upremas.F:75:72:

       75 |      &                             ierr)
   | 1
   Error: There is no specific subroutine for the generic
   'matmpiaijsetpreallocation' at (1)
   upremas.F:78:72:

       78 |      &                             ierr)
   | 1
   Error: There is no specific subroutine for the generic
   'matseqaijsetpreallocation' at (1)

--

On 3/23/25 1:09 AM, Jose E. Roman wrote:
"use petscmat" will use those definitions.

As I said, you probably have mismatching arguments. For instance
             call MatSeqAIJSetPreallocation(Mmat,
      &                             PETSC_NULL_INTEGER_ARRAY,mr(np(246)),
      &                             ierr)
The second argument is a PetscInt so PETSC_NULL_INTEGER_ARRAY is wrong, it 
should be PETSC_NULL_INTEGER.

Now the compiler will help you fix this kind of errors, which would go 
unnoticed before.

Jose



El 23 mar 2025, a las 9:02, Sanjay Govindjee<s...@berkeley.edu> escribió:

Jose,
   What module should I be using to load petscmat.h90?
-sanjay

--

On 3/23/25 12:54 AM, Jose E. Roman wrote:
Have a look at the list of changes - it is currently 
herehttps://petsc.org/main/changes/dev/ until the new version is released. See the last 
section "Fortran".

The functions ending in "F90" have been renamed, just remove the "F90" suffix.

Regarding the info-related errors, a workaround is to append %v, for instance
        mal = info(MAT_INFO_MALLOCS%v)
But Barry may want to provide a better fix for this.

Jose


El 23 mar 2025, a las 8:42, Jose E. Roman via 
petsc-users<petsc-users@mcs.anl.gov> escribió:

The Fortran interfaces for those functions are generated correctly, see 
$PETSC_ARCH/ftn/mat/petscmat.h90

For instance:

  interface MatMPIBAIJSetPreallocation
  subroutine MatMPIBAIJSetPreallocation(a,b,c,d,e,f, z)
  import tMat
  Mat :: a
  PetscInt :: b
  PetscInt :: c
  PetscInt :: d(*)
  PetscInt :: e
  PetscInt :: f(*)
  PetscErrorCode z
  end subroutine
  end interface

The compiler message is probably due to the type of an argument not matching 
the expected one. In particular, if you are passing NULL in one of the array 
arguments, you should use PETSC_NULL_INTEGER_ARRAY and not PETSC_NULL_INTEGER.

Jose


El 23 mar 2025, a las 8:25, Sanjay Govindjee via 
petsc-users<petsc-users@mcs.anl.gov> escribió:

Small update.  I managed to eliminate all the errors associated with 
PetscViewer and below (it had to do with the fact that I had not yet built a 
module that was needed).  The errors related to the preallocation routines 
still persists.
-sanjay

On 3/23/25 12:19 AM, Sanjay Govindjee wrote:
Hi Barry,
  I have moved to main and rebuilt the PETSc libraries etc.  Right now I am 
having trouble just getting my source code to compile.   Plenty of subroutines 
with PETSc calls compile but a few are throwing errors and killing my compile. 
I suspect there will be more but if I can figure these hopefully I can debug 
the ones that will follow.
-sanjay
Error: There is no specific subroutine for the generic 
'matmpibaijsetpreallocation' at (1)
upremas.F:68:72:

   68 |      &                               ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 
'matseqbaijsetpreallocation' at (1)
upremas.F:74:72:

   74 |      &                             ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 
'matmpiaijsetpreallocation' at (1)
upremas.F:77:72:

   77 |      &                             ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 
'matseqaijsetpreallocation' at (1)

parkv.F:58:25:

   58 |       PetscViewer    Y_view
      |                         1
Error: Type name 'tpetscviewer' at (1) is ambiguous
parkv.F:69:9:

   69 |       endif
      |         1
Error: Expecting END SUBROUTINE statement at (1)
parkv.F:72:9:

   72 |       endif
      |         1
Error: Expecting END SUBROUTINE statement at (1)
parkv.F:91:66:

   91 |         call PetscViewerASCIIOpen(PETSC_COMM_WORLD,"yvec.m",Y_view,
      |                                                                  1
Error: Symbol 'y_view' at (1) has no IMPLICIT type; did you mean 'yvec'?
parkv.F:65:72:

   65 |         call VecCreate        (PETSC_COMM_WORLD, xvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'veccreate' at (1)
parkv.F:67:72:

   67 |         call VecSetFromOptions(xvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecsetfromoptions' at 
(1)
parkv.F:68:72:

   68 |         call VecDuplicate     (xvec, yvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecduplicate' at (1)
parkv.F:71:72:

   71 |         call VecDuplicate     (xvec, yvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecduplicate' at (1)
parkv.F:85:72:

   85 |       call VecAssemblyBegin(xvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecassemblybegin' at (1)
parkv.F:86:72:

   86 |       call VecAssemblyEnd  (xvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecassemblyend' at (1)
parkv.F:88:72:

   88 |       call MatMult           (Kmat, xvec, yvec, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'matmult' at (1)
parkv.F:101:72:

  101 |       call VecGetOwnershipRange(yvec, starti, endi, ierr)
      |                                                                        1
Error: There is no specific subroutine for the generic 'vecgetownershiprange' 
at (1)


-


On 3/21/25 7:17 AM, Barry Smith wrote:
    I have just pushed a major update to the Fortran interface to the main 
PETSc git branch. Could you please try to work with main (to become release in 
a couple of weeks) with your Fortran code as we debug the problem? This will 
save you a lot of work and hopefully make the debugging more straightforward.

    You can send the same output with the debugger if it crashes in the main 
branch and I can try to track down what is going wrong.

  Barry




On Mar 21, 2025, at 12:37 AM, Sanjay Govindjee via 
petsc-users<petsc-users@mcs.anl.gov> wrote:

I am trying to upgrade my code to PETSc 3.22.4 (the code was last updated to 
3.19.4 or perhaps 3.18.1, I've lost track). I've been using this code with 
PETSc for over 20 years.

To get my code to compile and link during this update, I only need to make two 
changes; one was to use PetscViewerPushFormat instead of PetscViewerSetFormat 
and the other was to use PETSC_NULL_INTEGER_ARRAY in a spot or two.

When I run the code however, I am getting an error very early on during a call 
to MatCreate near the beginning of the code.  The screen output says:
[3]PETSC ERROR: matcreate_() at 
/Users/sg/petsc-3.22.4/gnug/src/mat/utils/ftn-auto/gcreatef.c:101 Cannot create 
PETSC_NULL_XXX object
[0]PETSC ERROR: matcreate_() at 
/Users/sg/petsc-3.22.4/gnug/src/mat/utils/ftn-auto/gcreatef.c:101 Cannot create 
PETSC_NULL_XXX object
[1]PETSC ERROR: matcreate_() at 
/Users/sg/petsc-3.22.4/gnug/src/mat/utils/ftn-auto/gcreatef.c:101 Cannot create 
PETSC_NULL_XXX object
[2]PETSC ERROR: matcreate_() at 
/Users/sg/petsc-3.22.4/gnug/src/mat/utils/ftn-auto/gcreatef.c:101 Cannot create 
PETSC_NULL_XXX object
I have a 4 processor run going.  I am running with -on_error_attach_debugger 
but the debugger is giving me cryptic (at least to me) output (the same for all 
4 processes modulo the PID).  Stack traces seem to be unavailable :(
lldb  -p 71963
(lldb) process attach --pid 71963
Process 71963 stopped
* thread #1, queue = 'com.apple.main-thread', stop reason = signal SIGSTOP
    frame #0: 0x00007fff69d92746 libsystem_kernel.dylib`__semwait_signal + 10
libsystem_kernel.dylib`__semwait_signal:
->  0x7fff69d92746 <+10>: jae    0x7fff69d92750            ; <+20>
    0x7fff69d92748 <+12>: movq   %rax, %rdi
    0x7fff69d9274b <+15>: jmp    0x7fff69d9121d            ; cerror
    0x7fff69d92750 <+20>: retq
Target 0: (feap) stopped.

Executable module set to "/Users/sg/Feap/ver87/parfeap/feap".
Architecture set to: x86_64h-apple-macosx-.
Does anyone have any hints as to what may be going on?  Note the program starts 
normally and i can do stuff with the interactive interface for the code -- even 
plotting the mesh etc. so I believe the input data has been read in correctly.  
The crash only occurs when I initiate the formation of the matrix.

I am attaching the 
/Users/sg/petsc-3.22.4/gnug/src/mat/utils/ftn-auto/gcreatef.c file in case that 
offers some insight.

Note, I have been
-sanjay
--

<gcreatef.c>
!$Id:$
      subroutine upremas(fl)

!      * * F E A P * * A Finite Element Analysis Program

!....  Copyright (c) 1984-2025: Regents of the University of California
!                               All rights reserved

!-----[--.----+----.----+----.-----------------------------------------]
!     Modification log                                Date (dd/mm/year)
!       Original version                                    01/11/2006
!       1. Added petscksp.h                                 01/05/2007
!       2. Change 'include/finclude' to 'finclude'          23/01/2009
!       3. Remove common 'pfeapa' (values in 'setups')      05/02/2009
!       4. Update petsc includes to v3.1                    20/07/2010
!       5. Change 'ndf' to 'nbsk' to match usolve.F         05/01/2013
!       6. Update for loss of VecValid, MatValid            05/01/2013
!       7. Update matrix creation calls                     06/01/2013
!       8. finclude -> petsc/finclude                       12/05/2016
!       9. Update for PETSc 3.8.3                           28/02/2018
!      10. Update pfeapc to a module                        28/03/2018
!      11. Remove unused 'chk'                              23/05/2019
!      12. Fixed NULL_INT to DEFAULT_INT for preallocation  04/06/2020
!      13. Update for PETSc 3.22.4                          20/03/2025
!-----[--.----+----.----+----.-----------------------------------------]
!     Purpose:  Mass interface for PETSc

!     Inputs:
!        fl(1) - Form Consistent mass if true
!        fl(2) - Form Lumped mass if true

!     Outputs:
!-----[--.----+----.----+----.-----------------------------------------]
#     include   <petsc/finclude/petscksp.h>
      use                       petscksp
      use                       petscmat   ! Newly added to try to fix compile 
errors
      use                       pfeapc
      implicit   none

#     include   "cdata.h"
#     include   "comblk.h"
#     include   "endata.h"
#     include   "iofile.h"
#     include   "pointer.h"
#     include   "sdata.h"
#     include   "setups.h"
#     include   "pfeapb.h"

      PetscErrorCode :: ierr

      logical        :: fl(2)

!     Preform setup
      if(fl(1)) then   ! Consistent mass allocate
        if(Mmat.eq.PETSC_NULL_MAT) then
          call MatCreate(PETSC_COMM_WORLD,Mmat,ierr)
          call MatSetSizes(Mmat,numpeq,numpeq,PETSC_DETERMINE, 
     &                        PETSC_DETERMINE,ierr)
          if(pfeap_bcin) call MatSetBlockSize(Mmat,nsbk,ierr)
          call MatSetFromOptions(Mmat, ierr)
          if(pfeap_blk) then
            call MatSetType(Mmat,MATBAIJ,ierr)
            call MatMPIBAIJSetPreallocation(Mmat,nsbk, 
     &                               PETSC_NULL_INTEGER,mr(np(246)), 
     &                               PETSC_NULL_INTEGER,mr(np(247)), 
     &                               ierr)
            call MatSeqBAIJSetPreallocation(Mmat,nsbk, 
     &                               PETSC_DEFAULT_INTEGER,mr(np(246)), 
     &                               ierr)
          else
            call MatSetType(Mmat,MATAIJ,ierr)
            call MatMPIAIJSetPreallocation(Mmat, 
     &                             PETSC_NULL_INTEGER,mr(np(246)),
     &                             PETSC_NULL_INTEGER,mr(np(247)),
     &                             ierr)
            call MatSeqAIJSetPreallocation(Mmat, 
     &                             PETSC_NULL_INTEGER,mr(np(246)),
     &                             ierr)
          endif
        endif
      elseif(fl(2)) then   ! Lumped mass allocate

        if(Mdiag.eq.PETSC_NULL_VEC) then
          call VecCreate        (PETSC_COMM_WORLD, Mdiag, ierr)
          call VecSetSizes      (Mdiag, numpeq, PETSC_DECIDE, ierr)
          call VecSetFromOptions(Mdiag, ierr)
        endif
      elseif(.not.fl(1) .and. .not.fl(2)) then
        write(iow,*) ' ERROR DID NOT ALLOCATE MASS MATRIX'
        write(ilg,*) ' ERROR DID NOT ALLOCATE MASS MATRIX'
        if(rank.eq.0) then
          write(*,*) ' ERROR DID NOT ALLOCATE MASS MATRIX'
        endif
        call plstop(.true.)
      endif
      if(ierr .ne. 0) then
        write(iow,*) 'Error on MatCreate'
        write(ilg,*) 'Error on MatCreate'
        if(rank.eq.0) then
          write(  *,*) 'Error on MatCreate'
        endif
        call plstop(.true.)
      endif

      if(fl(1)) then
        call MatZeroEntries   (Mmat,ierr)
      else
        call VecZeroEntries   (Mdiag,ierr)
      endif

!     Vector for matrix multiply
      if(xvec.eq.PETSC_NULL_VEC) then
        call VecCreate        (PETSC_COMM_WORLD, xvec, ierr)
        call VecSetSizes      (xvec, numpeq, PETSC_DECIDE, ierr)
        call VecSetFromOptions(xvec, ierr)
      endif

      end subroutine upremas

Reply via email to