[petsc-users] Problems with DMDAVecGetArrayF90 + Intel

Satish Balay balay at mcs.anl.gov
Tue May 15 11:33:21 CDT 2018


On Tue, 15 May 2018, Jed Brown wrote:

> Satish Balay <balay at mcs.anl.gov> writes:
> 
> > On Tue, 15 May 2018, Jed Brown wrote:
> >
> >> Wow, how did this ever work with other compilers?
> >> 
> >> To ensure everything gets fixed, I see we have
> >> 
> >>   #define F90Array3d void
> >> 
> >> and then arguments are defined as
> >> 
> >>   F90Array3d *ptr
> >> 
> >> We could make this type-safe by
> >> 
> >>   typedef struct { void *ptr; } F90Array3d;
> >
> > The following appears to be sufficient. [void* more appropriate than void - as these are pointers anyway?]
> 
> They're declared as pointers, so it's "pointer to void" versus "pointer
> to void*".

Your notation requires the following change

PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayf901_(DM *da,Vec *v,F90Array1d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
PETSC_EXTERN PetscErrorCode F90Array1dCreate(void*,MPI_Datatype,PetscInt,PetscInt,F90Array1d* PETSC_F90_2PTR_PROTO_NOVAR);

to:

PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayf901_(DM *da,Vec *v,F90Array1d a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
PETSC_EXTERN PetscErrorCode F90Array1dCreate(void*,MPI_Datatype,PetscInt,PetscInt,F90Array1d PETSC_F90_2PTR_PROTO_NOVAR);

[this doesn't look right.] So I guess 'pointer to address (void*)' is
the correct notation - as we pass these back to fortran calling
routine wrt dmdavecgetarrayf901_() and an output argument for
F90Array1dCreate().

Satish

> But the point of my suggestion was to get stronger type
> checking, including between arrays of different dimensions.
> 
> > -#define F90Array1d void
> > -#define F90Array2d void
> > -#define F90Array3d void
> > -#define F90Array4d void
> > +typedef void* F90Array1d;
> > +typedef void* F90Array2d;
> > +typedef void* F90Array3d;
> > +typedef void* F90Array4d;
> >
> >> 
> >> and use this for arguments:
> >> 
> >>   F90Array3d a
> >
> > Don't need this change..
> >
> > Attaching my current patch..
> >
> > Satish
> > diff --git a/include/petsc/private/f90impl.h b/include/petsc/private/f90impl.h
> > index a35efb76bd..4f26c8ffea 100644
> > --- a/include/petsc/private/f90impl.h
> > +++ b/include/petsc/private/f90impl.h
> > @@ -15,11 +15,10 @@
> >  #endif
> >  
> >  #if defined(PETSC_USING_F90)
> > -
> > -#define F90Array1d void
> > -#define F90Array2d void
> > -#define F90Array3d void
> > -#define F90Array4d void
> > +typedef void* F90Array1d;
> > +typedef void* F90Array2d;
> > +typedef void* F90Array3d;
> > +typedef void* F90Array4d;
> >  
> >  PETSC_EXTERN PetscErrorCode F90Array1dCreate(void*,MPI_Datatype,PetscInt,PetscInt,F90Array1d* PETSC_F90_2PTR_PROTO_NOVAR);
> >  PETSC_EXTERN PetscErrorCode F90Array1dAccess(F90Array1d*,MPI_Datatype,void** PETSC_F90_2PTR_PROTO_NOVAR);
> > diff --git a/src/dm/impls/composite/f90-custom/zfddaf90.c b/src/dm/impls/composite/f90-custom/zfddaf90.c
> > index abc0a27a01..2633c9e86a 100644
> > --- a/src/dm/impls/composite/f90-custom/zfddaf90.c
> > +++ b/src/dm/impls/composite/f90-custom/zfddaf90.c
> > @@ -24,8 +24,8 @@ PETSC_EXTERN void PETSC_STDCALL dmcompositegetaccessvpvp_(DM *dm,Vec *v,Vec *v1,
> >  PETSC_EXTERN void PETSC_STDCALL dmcompositerestoreaccessvpvp_(DM *dm,Vec *v,Vec *v1,F90Array1d *p1,Vec *v2,F90Array1d *p2,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd1) PETSC_F90_2PTR_PROTO(ptrd2))
> >  {
> >    *ierr = DMCompositeRestoreAccess(*dm,*v,v1,0,v2,0);
> > -  *ierr = F90Array1dDestroy(&p1,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd1));
> > -  *ierr = F90Array1dDestroy(&p2,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd2));
> > +  *ierr = F90Array1dDestroy(p1,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd1));
> > +  *ierr = F90Array1dDestroy(p2,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd2));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmcompositegetentriesarray_(DM *dm, DM *dmarray, PetscErrorCode *ierr)
> > diff --git a/src/dm/impls/da/f90-custom/zda1f90.c b/src/dm/impls/da/f90-custom/zda1f90.c
> > index 082027725f..41cc58534f 100644
> > --- a/src/dm/impls/da/f90-custom/zda1f90.c
> > +++ b/src/dm/impls/da/f90-custom/zda1f90.c
> > @@ -74,7 +74,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayf901_(DM *da,Vec *v,F90Array1
> >    PetscScalar *fa;
> >    *ierr = F90Array1dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArray(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array1dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array1dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayf902_(DM *da,Vec *v,F90Array2d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -113,7 +113,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayf902_(DM *da,Vec *v,F90Array2
> >    PetscScalar *fa;
> >    *ierr = F90Array2dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArray(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array2dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array2dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayf903_(DM *da,Vec *v,F90Array3d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -154,7 +154,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayf903_(DM *da,Vec *v,F90Array3
> >    PetscScalar *fa;
> >    *ierr = F90Array3dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArray(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array3dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array3dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayf904_(DM *da,Vec *v,F90Array4d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -190,7 +190,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayf904_(DM *da,Vec *v,F90Array4
> >    */
> >    *ierr = F90Array4dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArray(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array4dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array4dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayreadf901_(DM *da,Vec *v,F90Array1d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -223,7 +223,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayreadf901_(DM *da,Vec *v,F90Ar
> >    const PetscScalar *fa;
> >    *ierr = F90Array1dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArrayRead(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array1dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array1dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayreadf902_(DM *da,Vec *v,F90Array2d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -262,7 +262,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayreadf902_(DM *da,Vec *v,F90Ar
> >    const PetscScalar *fa;
> >    *ierr = F90Array2dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArrayRead(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array2dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array2dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayreadf903_(DM *da,Vec *v,F90Array3d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -303,7 +303,7 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayreadf903_(DM *da,Vec *v,F90Ar
> >    const PetscScalar *fa;
> >    *ierr = F90Array3dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArrayRead(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array3dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array3dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> >  
> >  PETSC_EXTERN void PETSC_STDCALL dmdavecgetarrayreadf904_(DM *da,Vec *v,F90Array4d *a,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptrd))
> > @@ -339,5 +339,5 @@ PETSC_EXTERN void PETSC_STDCALL dmdavecrestorearrayreadf904_(DM *da,Vec *v,F90Ar
> >    */
> >    *ierr = F90Array4dAccess(a,MPIU_SCALAR,(void**)&fa PETSC_F90_2PTR_PARAM(ptrd));
> >    *ierr = VecRestoreArrayRead(*v,&fa);if (*ierr) return;
> > -  *ierr = F90Array4dDestroy(&a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> > +  *ierr = F90Array4dDestroy(a,MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
> >  }
> 



More information about the petsc-users mailing list