xref: /petsc/src/ksp/pc/impls/asm/ftn-custom/zasmf.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
1*af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
2c6db04a5SJed Brown #include <petscksp.h>
3e54e4138SSatish Balay 
4e54e4138SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
5e54e4138SSatish Balay #define pcasmgetsubksp_            PCASMGETSUBKSP
6e54e4138SSatish Balay #define pcasmsetlocalsubdomains_   PCASMSETLOCALSUBDOMAINS
7e54e4138SSatish Balay #define pcasmsetglobalsubdomains_  PCASMSETGLOBALSUBDOMAINS
8e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_  PCASMGETLOCALSUBMATRICES
9e54e4138SSatish Balay #define pcasmgetlocalsubdomains_   PCASMGETLOCALSUBDOMAINS
10ef356e52SBarry Smith #define pcasmcreatesubdomains_     PCASMCREATESUBDOMAINS
11ef356e52SBarry Smith #define pcasmdestroysubdomains_    PCASMDESTROYSUBDOMAINS
12e54e4138SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
13e54e4138SSatish Balay #define pcasmgetsubksp_            pcasmgetsubksp
14e54e4138SSatish Balay #define pcasmsetlocalsubdomains_   pcasmsetlocalsubdomains
15e54e4138SSatish Balay #define pcasmsetglobalsubdomains_  pcasmsetglobalsubdomains
16e54e4138SSatish Balay #define pcasmgetlocalsubmatrices_  pcasmgetlocalsubmatrices
17e54e4138SSatish Balay #define pcasmgetlocalsubdomains_   pcasmgetlocalsubdomains
18ef356e52SBarry Smith #define pcasmcreatesubdomains_     pcasmcreatesubdomains
19ef356e52SBarry Smith #define pcasmdestroysubdomains_    pcasmdestroysubdomains
20e54e4138SSatish Balay #endif
21e54e4138SSatish Balay 
228cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmcreatesubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
23ef356e52SBarry Smith {
24ef356e52SBarry Smith   PetscInt i;
25ef356e52SBarry Smith   IS       *insubs;
26ef356e52SBarry Smith 
27ef356e52SBarry Smith   *ierr = PCASMCreateSubdomains(*mat,*n,&insubs);if (*ierr) return;
282fa5cd67SKarl Rupp   for (i=0; i<*n; i++) subs[i] = insubs[i];
29ef356e52SBarry Smith   *ierr = PetscFree(insubs);
30ef356e52SBarry Smith }
31ef356e52SBarry Smith 
32ef356e52SBarry Smith 
338cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmdestroysubdomains_(Mat *mat,PetscInt *n,IS *subs,PetscErrorCode *ierr)
34ef356e52SBarry Smith {
35ef356e52SBarry Smith   PetscInt i;
36ef356e52SBarry Smith 
37ef356e52SBarry Smith   for (i=0; i<*n; i++) {
38fcfd50ebSBarry Smith     *ierr = ISDestroy(&subs[i]);if (*ierr) return;
39ef356e52SBarry Smith   }
40ef356e52SBarry Smith }
41ef356e52SBarry Smith 
428cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmgetsubksp_(PC *pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp,PetscErrorCode *ierr)
43e54e4138SSatish Balay {
44e54e4138SSatish Balay   KSP      *tksp;
45e54e4138SSatish Balay   PetscInt i,nloc;
46e54e4138SSatish Balay   CHKFORTRANNULLINTEGER(n_local);
47e54e4138SSatish Balay   CHKFORTRANNULLINTEGER(first_local);
48d29017ddSJed Brown   CHKFORTRANNULLOBJECT(ksp);
49e54e4138SSatish Balay   *ierr = PCASMGetSubKSP(*pc,&nloc,first_local,&tksp);
50e54e4138SSatish Balay   if (n_local) *n_local = nloc;
51d29017ddSJed Brown   if (ksp) {
522fa5cd67SKarl Rupp     for (i=0; i<nloc; i++) ksp[i] = tksp[i];
53e54e4138SSatish Balay   }
54d29017ddSJed Brown }
55e54e4138SSatish Balay 
568cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmsetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
57e54e4138SSatish Balay {
58e54e4138SSatish Balay   CHKFORTRANNULLOBJECT(is);
592b691e39SMatthew Knepley   CHKFORTRANNULLOBJECT(is_local);
602b691e39SMatthew Knepley   *ierr = PCASMSetLocalSubdomains(*pc,*n,is,is_local);
61e54e4138SSatish Balay }
62e54e4138SSatish Balay 
638cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmsettotalsubdomains_(PC *pc,PetscInt *N,IS *is,IS *is_local, PetscErrorCode *ierr)
64e54e4138SSatish Balay {
65e54e4138SSatish Balay   CHKFORTRANNULLOBJECT(is);
662b691e39SMatthew Knepley   CHKFORTRANNULLOBJECT(is_local);
672b691e39SMatthew Knepley   *ierr = PCASMSetTotalSubdomains(*pc,*N,is,is_local);
68e54e4138SSatish Balay }
69e54e4138SSatish Balay 
708cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmgetlocalsubmatrices_(PC *pc,PetscInt *n,Mat *mat, PetscErrorCode *ierr)
71e54e4138SSatish Balay {
72e54e4138SSatish Balay   PetscInt nloc,i;
73e54e4138SSatish Balay   Mat      *tmat;
74e54e4138SSatish Balay   CHKFORTRANNULLOBJECT(mat);
75e54e4138SSatish Balay   CHKFORTRANNULLINTEGER(n);
76e54e4138SSatish Balay   *ierr = PCASMGetLocalSubmatrices(*pc,&nloc,&tmat);
77e54e4138SSatish Balay   if (n) *n = nloc;
78e54e4138SSatish Balay   if (mat) {
792fa5cd67SKarl Rupp     for (i=0; i<nloc; i++) mat[i] = tmat[i];
80e54e4138SSatish Balay   }
81e54e4138SSatish Balay }
828cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL pcasmgetlocalsubdomains_(PC *pc,PetscInt *n,IS *is,IS *is_local, PetscErrorCode *ierr)
83e54e4138SSatish Balay {
84e54e4138SSatish Balay   PetscInt nloc,i;
852b691e39SMatthew Knepley   IS       *tis, *tis_local;
86e54e4138SSatish Balay   CHKFORTRANNULLOBJECT(is);
872b691e39SMatthew Knepley   CHKFORTRANNULLOBJECT(is_local);
88e54e4138SSatish Balay   CHKFORTRANNULLINTEGER(n);
892b691e39SMatthew Knepley   *ierr = PCASMGetLocalSubdomains(*pc,&nloc,&tis,&tis_local);
90e54e4138SSatish Balay   if (n) *n = nloc;
91e54e4138SSatish Balay   if (is) {
922fa5cd67SKarl Rupp     for (i=0; i<nloc; i++) is[i] = tis[i];
93e54e4138SSatish Balay   }
942b691e39SMatthew Knepley   if (is_local && tis_local) {
952fa5cd67SKarl Rupp     for (i=0; i<nloc; i++) is[i] = tis_local[i];
962b691e39SMatthew Knepley   }
97e54e4138SSatish Balay }
98e54e4138SSatish Balay 
99