xref: /petsc/src/ksp/pc/impls/fieldsplit/ftn-custom/zfieldsplitf.c (revision c3499a38069e8bafdf09d701215e611c558ad9a9)
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
9*c3499a38SBarry 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
15*c3499a38SBarry 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);
23285fb4e2SStefano Zampini   *ierr = PCFieldSplitSchurGetSubKSP(*pc,&nloc,&tksp); if (*ierr) return;
24285fb4e2SStefano Zampini   if (n_local) *n_local = nloc;
25285fb4e2SStefano Zampini   CHKFORTRANNULLOBJECT(ksp);
26285fb4e2SStefano Zampini   if (ksp) {
27285fb4e2SStefano Zampini     for (i=0; i<nloc; i++) ksp[i] = tksp[i];
28285fb4e2SStefano Zampini   }
29285fb4e2SStefano Zampini   *ierr = PetscFree(tksp);
30285fb4e2SStefano Zampini }
31285fb4e2SStefano Zampini 
3219caf8f3SSatish Balay PETSC_EXTERN void pcfieldsplitgetsubksp_(PC *pc,PetscInt *n_local,KSP *ksp,PetscErrorCode *ierr)
331193e19dSBarry Smith {
341193e19dSBarry Smith   KSP      *tksp;
351193e19dSBarry Smith   PetscInt i,nloc;
361193e19dSBarry Smith   CHKFORTRANNULLINTEGER(n_local);
371193e19dSBarry Smith   *ierr = PCFieldSplitGetSubKSP(*pc,&nloc,&tksp); if (*ierr) return;
381193e19dSBarry Smith   if (n_local) *n_local = nloc;
391193e19dSBarry Smith   CHKFORTRANNULLOBJECT(ksp);
401193e19dSBarry Smith   if (ksp) {
412fa5cd67SKarl Rupp     for (i=0; i<nloc; i++) ksp[i] = tksp[i];
421193e19dSBarry Smith   }
43e06cc25fSBarry Smith   *ierr = PetscFree(tksp);
441193e19dSBarry Smith }
451193e19dSBarry Smith 
4619caf8f3SSatish Balay PETSC_EXTERN void  pcfieldsplitsetis_(PC *pc, char* splitname,IS *is, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
47218d4943SBarry Smith {
48218d4943SBarry Smith   char *t;
49218d4943SBarry Smith   FIXCHAR(splitname,len,t);
50d49bb8f9SBarry Smith   *ierr = PCFieldSplitSetIS(*pc,t,*is);if (*ierr) return;
51218d4943SBarry Smith   FREECHAR(splitname,t);
52218d4943SBarry Smith }
53218d4943SBarry Smith 
5419caf8f3SSatish Balay PETSC_EXTERN void  pcfieldsplitgetis_(PC *pc, char* splitname,IS *is, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
55d6f11fddSBarry Smith {
56d6f11fddSBarry Smith   char *t;
57d6f11fddSBarry Smith   FIXCHAR(splitname,len,t);
58d6f11fddSBarry Smith   *ierr = PCFieldSplitGetIS(*pc,t,is);if (*ierr) return;
59d6f11fddSBarry Smith   FREECHAR(splitname,t);
60d6f11fddSBarry Smith }
61d6f11fddSBarry Smith 
62*c3499a38SBarry 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)
63*c3499a38SBarry Smith {
64*c3499a38SBarry Smith   char *t;
65*c3499a38SBarry Smith   FIXCHAR(splitname,len,t);
66*c3499a38SBarry Smith   *ierr = PCFieldSplitSetFields(*pc, t, *n, fields, fields_col);if (*ierr) return;
67*c3499a38SBarry Smith   FREECHAR(splitname,t);
68*c3499a38SBarry Smith }
69