xref: /petsc/src/mat/impls/nest/ftn-custom/zmatnestf.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
245c38901SJed Brown #include <petscmat.h>
345c38901SJed Brown 
445c38901SJed Brown #if defined(PETSC_HAVE_FORTRAN_CAPS)
545c38901SJed Brown #define matcreatenest_                   MATCREATENEST
63a4d7b9aSSatish Balay #define matnestgetiss_                   MATNESTGETISS
7ffa9b3b1SVincent Le Chenadec #define matnestgetsubmats_               MATNESTGETSUBMATS
845c38901SJed Brown #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
945c38901SJed Brown #define matcreatenest_                   matcreatenest
103a4d7b9aSSatish Balay #define matnestgetiss_                   matnestgetiss
11ffa9b3b1SVincent Le Chenadec #define matnestgetsubmats_               matnestgetsubmats
1245c38901SJed Brown #endif
1345c38901SJed Brown 
148cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL matcreatenest_(MPI_Fint *comm,PetscInt *nr,IS is_row[],PetscInt *nc,IS is_col[],Mat a[],Mat *B,int *ierr)
1545c38901SJed Brown {
16*2f6eced2SAlex Fikl   Mat *m,*tmp;
17*2f6eced2SAlex Fikl   PetscInt i;
18*2f6eced2SAlex Fikl 
1945c38901SJed Brown   CHKFORTRANNULLOBJECT(is_row);
2045c38901SJed Brown   CHKFORTRANNULLOBJECT(is_col);
21*2f6eced2SAlex Fikl 
22*2f6eced2SAlex Fikl   *ierr = PetscMalloc1((*nr)*(*nc), &m); if(*ierr) return;
23*2f6eced2SAlex Fikl   for(i=0; i<(*nr)*(*nc); i++) {
24*2f6eced2SAlex Fikl     tmp = &(a[i]);
25*2f6eced2SAlex Fikl     CHKFORTRANNULLOBJECT(tmp);
26*2f6eced2SAlex Fikl     m[i] = (tmp == NULL ? NULL : a[i]);
27*2f6eced2SAlex Fikl   }
28*2f6eced2SAlex Fikl   *ierr = MatCreateNest(MPI_Comm_f2c(*comm),*nr,is_row,*nc,is_col,m,B); if (*ierr) return;
29*2f6eced2SAlex Fikl   *ierr = PetscFree(m);
3045c38901SJed Brown }
313a4d7b9aSSatish Balay 
323a4d7b9aSSatish Balay PETSC_EXTERN void PETSC_STDCALL  matnestgetiss_(Mat *A,IS rows[],IS cols[], int *ierr )
333a4d7b9aSSatish Balay {
343a4d7b9aSSatish Balay   CHKFORTRANNULLOBJECT(rows);
353a4d7b9aSSatish Balay   CHKFORTRANNULLOBJECT(cols);
363a4d7b9aSSatish Balay   *ierr = MatNestGetISs(*A,rows,cols);
373a4d7b9aSSatish Balay }
38ffa9b3b1SVincent Le Chenadec 
39ffa9b3b1SVincent Le Chenadec PETSC_EXTERN void PETSC_STDCALL matnestgetsubmats_(Mat *A,PetscInt *M,PetscInt *N,Mat *sub,int *ierr)
40ffa9b3b1SVincent Le Chenadec {
41351962e3SVincent Le Chenadec   PetscInt i,j,m,n;
42ffa9b3b1SVincent Le Chenadec   Mat **mat;
43351962e3SVincent Le Chenadec 
44351962e3SVincent Le Chenadec   CHKFORTRANNULLINTEGER(M);
45351962e3SVincent Le Chenadec   CHKFORTRANNULLINTEGER(N);
46351962e3SVincent Le Chenadec   CHKFORTRANNULLOBJECT(sub);
47351962e3SVincent Le Chenadec 
48351962e3SVincent Le Chenadec   *ierr = MatNestGetSubMats(*A,&m,&n,&mat);
49351962e3SVincent Le Chenadec 
50351962e3SVincent Le Chenadec   if (M) {
51351962e3SVincent Le Chenadec     *M = m;
52351962e3SVincent Le Chenadec   }
53351962e3SVincent Le Chenadec   if (N) {
54351962e3SVincent Le Chenadec     *N = n;
55351962e3SVincent Le Chenadec   }
56351962e3SVincent Le Chenadec   if (sub) {
57351962e3SVincent Le Chenadec     for (i=0; i<m; i++) {
58351962e3SVincent Le Chenadec       for (j=0; j<n; j++) {
59*2f6eced2SAlex Fikl         if (mat[i][j]) {
60351962e3SVincent Le Chenadec           sub[j + n * i] = mat[i][j];
61*2f6eced2SAlex Fikl         } else {
62*2f6eced2SAlex Fikl           sub[j + n * i] = (void *)-1;
63*2f6eced2SAlex Fikl         }
64351962e3SVincent Le Chenadec       }
65ffa9b3b1SVincent Le Chenadec     }
66ffa9b3b1SVincent Le Chenadec   }
67ffa9b3b1SVincent Le Chenadec }
68