xref: /petsc/src/ksp/pc/impls/fieldsplit/ftn-custom/zfieldsplitf.c (revision 5975b3b6e3931510e2a64a701673cbe1930c6f42)
1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
21193e19dSBarry Smith #include <petscksp.h>
31193e19dSBarry Smith 
41193e19dSBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
51193e19dSBarry Smith   #define pcfieldsplitgetsubksp_      PCFIELDSPLITGETSUBKSP
6285fb4e2SStefano Zampini   #define pcfieldsplitschurgetsubksp_ PCFIELDSPLITSCHURGETSUBKSP
7218d4943SBarry Smith   #define pcfieldsplitsetis_          PCFIELDSPLITSETIS
8d6f11fddSBarry Smith   #define pcfieldsplitgetis_          PCFIELDSPLITGETIS
9c3499a38SBarry Smith   #define pcfieldsplitsetfields_      PCFIELDSPLITSETFIELDS
101193e19dSBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
111193e19dSBarry Smith   #define pcfieldsplitgetsubksp_      pcfieldsplitgetsubksp
12285fb4e2SStefano Zampini   #define pcfieldsplitschurgetsubksp_ pcfieldsplitschurgetsubksp
13218d4943SBarry Smith   #define pcfieldsplitsetis_          pcfieldsplitsetis
14d6f11fddSBarry Smith   #define pcfieldsplitgetis_          pcfieldsplitgetis
15c3499a38SBarry Smith   #define pcfieldsplitsetfields_      pcfieldsplitsetfields
161193e19dSBarry Smith #endif
171193e19dSBarry Smith 
1819caf8f3SSatish Balay PETSC_EXTERN void pcfieldsplitschurgetsubksp_(PC *pc, PetscInt *n_local, KSP *ksp, PetscErrorCode *ierr)
19285fb4e2SStefano Zampini {
20285fb4e2SStefano Zampini   KSP     *tksp;
21285fb4e2SStefano Zampini   PetscInt i, nloc;
22285fb4e2SStefano Zampini   CHKFORTRANNULLINTEGER(n_local);
23*5975b3b6SBarry Smith   *ierr = PCFieldSplitSchurGetSubKSP(*pc, &nloc, &tksp);
24*5975b3b6SBarry Smith   if (*ierr) return;
25285fb4e2SStefano Zampini   if (n_local) *n_local = nloc;
26285fb4e2SStefano Zampini   CHKFORTRANNULLOBJECT(ksp);
27285fb4e2SStefano Zampini   if (ksp) {
28285fb4e2SStefano Zampini     for (i = 0; i < nloc; i++) ksp[i] = tksp[i];
29285fb4e2SStefano Zampini   }
30285fb4e2SStefano Zampini   *ierr = PetscFree(tksp);
31285fb4e2SStefano Zampini }
32285fb4e2SStefano Zampini 
3319caf8f3SSatish Balay PETSC_EXTERN void pcfieldsplitgetsubksp_(PC *pc, PetscInt *n_local, KSP *ksp, PetscErrorCode *ierr)
341193e19dSBarry Smith {
351193e19dSBarry Smith   KSP     *tksp;
361193e19dSBarry Smith   PetscInt i, nloc;
371193e19dSBarry Smith   CHKFORTRANNULLINTEGER(n_local);
38*5975b3b6SBarry Smith   *ierr = PCFieldSplitGetSubKSP(*pc, &nloc, &tksp);
39*5975b3b6SBarry Smith   if (*ierr) return;
401193e19dSBarry Smith   if (n_local) *n_local = nloc;
411193e19dSBarry Smith   CHKFORTRANNULLOBJECT(ksp);
421193e19dSBarry Smith   if (ksp) {
432fa5cd67SKarl Rupp     for (i = 0; i < nloc; i++) ksp[i] = tksp[i];
441193e19dSBarry Smith   }
45e06cc25fSBarry Smith   *ierr = PetscFree(tksp);
461193e19dSBarry Smith }
471193e19dSBarry Smith 
4819caf8f3SSatish Balay PETSC_EXTERN void pcfieldsplitsetis_(PC *pc, char *splitname, IS *is, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
49218d4943SBarry Smith {
50218d4943SBarry Smith   char *t;
51218d4943SBarry Smith   FIXCHAR(splitname, len, t);
52*5975b3b6SBarry Smith   *ierr = PCFieldSplitSetIS(*pc, t, *is);
53*5975b3b6SBarry Smith   if (*ierr) return;
54218d4943SBarry Smith   FREECHAR(splitname, t);
55218d4943SBarry Smith }
56218d4943SBarry Smith 
5719caf8f3SSatish Balay PETSC_EXTERN void pcfieldsplitgetis_(PC *pc, char *splitname, IS *is, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
58d6f11fddSBarry Smith {
59d6f11fddSBarry Smith   char *t;
60d6f11fddSBarry Smith   FIXCHAR(splitname, len, t);
61*5975b3b6SBarry Smith   *ierr = PCFieldSplitGetIS(*pc, t, is);
62*5975b3b6SBarry Smith   if (*ierr) return;
63d6f11fddSBarry Smith   FREECHAR(splitname, t);
64d6f11fddSBarry Smith }
65d6f11fddSBarry Smith 
66c3499a38SBarry Smith PETSC_EXTERN void pcfieldsplitsetfields_(PC *pc, char *splitname, PetscInt *n, const PetscInt *fields, const PetscInt *fields_col, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
67c3499a38SBarry Smith {
68c3499a38SBarry Smith   char *t;
69c3499a38SBarry Smith   FIXCHAR(splitname, len, t);
70*5975b3b6SBarry Smith   *ierr = PCFieldSplitSetFields(*pc, t, *n, fields, fields_col);
71*5975b3b6SBarry Smith   if (*ierr) return;
72c3499a38SBarry Smith   FREECHAR(splitname, t);
73c3499a38SBarry Smith }
74