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