xref: /petsc/src/vec/is/utils/ftn-custom/zvsectionisf.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 #include <petsc/private/fortranimpl.h>
2 #include <petscis.h>
3 #include <petscsection.h>
4 #include <petscviewer.h>
5 
6 #if defined(PETSC_HAVE_FORTRAN_CAPS)
7 #define petscsectiongetpointsyms_          PETSCSECTIONGETPOINTSYMS
8 #define petscsectionrestorepointsyms_      PETSCSECTIONRESTOREPOINTSYMS
9 #define petscsectiongetfieldpointsyms_     PETSCSECTIONGETFIELDPOINTSYMS
10 #define petscsectionrestorefieldpointsyms_ PETSCSECTIONRESTOREFIELDPOINTSYMS
11 #define petscsectionview_                  PETSCSECTIONVIEW
12 #define petscsectiongetfieldname_          PETSCSECTIONGETFIELDNAME
13 #define petscsectionsetfieldname_          PETSCSECTIONSETFIELDNAME
14 #define petscsfdistributesection_          PETSCSFDISTRIBUTESECTION
15 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
16 #define petscsectiongetpointsyms_          petscsectiongetpointsyms
17 #define petscsectionrestorepointsyms_      petscsectionrestorepointsyms
18 #define petscsectiongetfieldpointsyms_     petscsectiongetfieldpointsyms
19 #define petscsectionrestorefieldpointsyms_ petscsectionrestorefieldpointsyms
20 #define petscsectionview_                  petscsectionview
21 #define petscsectiongetfieldname_          petscsectiongetfieldname
22 #define petscsectionsetfieldname_          petscsectionsetfieldname
23 #define petscsfdistributesection_          petscsfdistributesection
24 #endif
25 
26 PETSC_EXTERN void  petscsectiongetpointsyms_(PetscSection section,PetscInt *numPoints, PetscInt *points, PetscInt ***perms, PetscScalar ***rots, int *__ierr ){
27 *__ierr = PetscSectionGetPointSyms(section,*numPoints,points,(const PetscInt ***)perms,(const PetscScalar ***)rots);
28 }
29 PETSC_EXTERN void  petscsectionrestorepointsyms_(PetscSection section,PetscInt *numPoints, PetscInt *points, PetscInt ***perms, PetscScalar ***rots, int *__ierr ){
30 *__ierr = PetscSectionRestorePointSyms(section,*numPoints,points,(const PetscInt ***)perms,(const PetscScalar ***)rots);
31 }
32 PETSC_EXTERN void  petscsectiongetfieldpointsyms_(PetscSection section,PetscInt *field,PetscInt *numPoints, PetscInt *points, PetscInt ***perms, PetscScalar ***rots, int *__ierr ){
33 *__ierr = PetscSectionGetFieldPointSyms(section,*field,*numPoints,points,(const PetscInt ***)perms,(const PetscScalar ***)rots);
34 }
35 PETSC_EXTERN void  petscsectionrestorefieldpointsyms_(PetscSection section,PetscInt *field,PetscInt *numPoints, PetscInt *points, PetscInt ***perms, PetscScalar ***rots, int *__ierr ){
36 *__ierr = PetscSectionRestoreFieldPointSyms(section,*field,*numPoints,points,(const PetscInt ***)perms,(const PetscScalar ***)rots);
37 }
38 
39 PETSC_EXTERN void petscsectionview_(PetscSection *s, PetscViewer *vin, PetscErrorCode *ierr)
40 {
41   PetscViewer v;
42 
43   PetscPatchDefaultViewers_Fortran(vin, v);
44   *ierr = PetscSectionView(*s, v);
45 }
46 
47 PETSC_EXTERN void petscsectiongetfieldname_(PetscSection *s, PetscInt *field, char* name, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
48 {
49   const char *fname;
50 
51   *ierr = PetscSectionGetFieldName(*s, *field, &fname);if (*ierr) return;
52   *ierr = PetscStrncpy(name, fname, len);
53   FIXRETURNCHAR(PETSC_TRUE,name,len);
54 }
55 
56 PETSC_EXTERN void petscsectionsetfieldname_(PetscSection *s, PetscInt *field, char* name, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
57 {
58   char *f;
59 
60   FIXCHAR(name, len, f);
61   *ierr = PetscSectionSetFieldName(*s, *field, f);if (*ierr) return;
62   FREECHAR(name, f);
63 }
64 
65 PETSC_EXTERN void  petscsfdistributesection_(PetscSF sf,PetscSection rootSection,PetscInt **remoteOffsets,PetscSection leafSection, int *__ierr ){
66   if (remoteOffsets != PETSC_NULL_INTEGER_Fortran) {
67     PetscError(PETSC_COMM_SELF, __LINE__, "PetscSFDistributeSection_Fortran", __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL,
68                "The remoteOffsets argument must be PETSC_NULL_INTEGER in Fortran");
69     *__ierr = 1;
70     return;
71   }
72   *__ierr = PetscSFDistributeSection(sf,rootSection,NULL,leafSection);
73 }
74