xref: /petsc/src/mat/graphops/partition/ftn-custom/zpartitionf.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 #include <../src/mat/impls/adj/mpi/mpiadj.h>
2 #include <petsc/private/fortranimpl.h>
3 #include <petscmat.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6   #define matpartitioningsetvertexweights_ MATPARTITIONINGSETVERTEXWEIGHTS
7   #define matpartitioningview_             MATPARTITIONINGVIEW
8   #define matpartitioningsettype_          MATPARTITIONINGSETTYPE
9   #define matpartitioningviewfromoptions_  MATPARTITIONINGVIEWFROMOPTIONS
10 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
11   #define matpartitioningsetvertexweights_ matpartitioningsetvertexweights
12   #define matpartitioningview_             matpartitioningview
13   #define matpartitioningsettype_          matpartitioningsettype
14   #define matpartitioningviewfromoptions_  matpartitioningviewfromoptions
15 #endif
16 
17 PETSC_EXTERN void matpartitioningsetvertexweights_(MatPartitioning *part, const PetscInt weights[], PetscErrorCode *ierr)
18 {
19   PetscInt  len;
20   PetscInt *array;
21   *ierr = MatGetLocalSize((*part)->adj, &len, NULL);
22   if (*ierr) return;
23   *ierr = PetscMalloc1(len, &array);
24   if (*ierr) return;
25   *ierr = PetscArraycpy(array, weights, len);
26   if (*ierr) return;
27   *ierr = MatPartitioningSetVertexWeights(*part, array);
28 }
29 PETSC_EXTERN void matpartitioningview_(MatPartitioning *part, PetscViewer *viewer, PetscErrorCode *ierr)
30 {
31   PetscViewer v;
32   PetscPatchDefaultViewers_Fortran(viewer, v);
33   *ierr = MatPartitioningView(*part, v);
34 }
35 
36 PETSC_EXTERN void matpartitioningsettype_(MatPartitioning *part, char *type, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
37 {
38   char *t;
39   FIXCHAR(type, len, t);
40   *ierr = MatPartitioningSetType(*part, t);
41   if (*ierr) return;
42   FREECHAR(type, t);
43 }
44 PETSC_EXTERN void matpartitioningviewfromoptions_(MatPartitioning *ao, PetscObject obj, char *type, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
45 {
46   char *t;
47 
48   FIXCHAR(type, len, t);
49   CHKFORTRANNULLOBJECT(obj);
50   *ierr = MatPartitioningViewFromOptions(*ao, obj, t);
51   if (*ierr) return;
52   FREECHAR(type, t);
53 }
54