Actual source code: aijmatlab.c

  1: /*
  2:         Provides an interface for the MATLAB engine sparse solver

  4: */
  5: #include <../src/mat/impls/aij/seq/aij.h>
  6: #include <petscmatlab.h>
  7: #include <engine.h> /* MATLAB include file */
  8: #include <mex.h>    /* MATLAB include file */

 10: static mxArray *MatSeqAIJToMatlab(Mat B)
 11: {
 12:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)B->data;
 13:   mwIndex    *ii, *jj;
 14:   mxArray    *mat;
 15:   PetscInt    i;

 17:   mat = mxCreateSparse(B->cmap->n, B->rmap->n, aij->nz, mxREAL);
 18:   if (PetscArraycpy(mxGetPr(mat), aij->a, aij->nz)) return NULL;
 19:   /* MATLAB stores by column, not row so we pass in the transpose of the matrix */
 20:   jj = mxGetIr(mat);
 21:   for (i = 0; i < aij->nz; i++) jj[i] = aij->j[i];
 22:   ii = mxGetJc(mat);
 23:   for (i = 0; i < B->rmap->n + 1; i++) ii[i] = aij->i[i];
 24:   return mat;
 25: }

 27: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject obj, void *mengine)
 28: {
 29:   mxArray *mat;

 31:   PetscFunctionBegin;
 32:   mat = MatSeqAIJToMatlab((Mat)obj);
 33:   PetscCheck(mat, PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot create MATLAB matrix");
 34:   PetscCall(PetscObjectName(obj));
 35:   engPutVariable((Engine *)mengine, obj->name, mat);
 36:   PetscFunctionReturn(PETSC_SUCCESS);
 37: }

 39: static PetscErrorCode MatSeqAIJFromMatlab(mxArray *mmat, Mat mat)
 40: {
 41:   PetscInt    nz, n, m, *i, *j, k;
 42:   mwIndex     nnz, nn, nm, *ii, *jj;
 43:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;

 45:   PetscFunctionBegin;
 46:   nn  = mxGetN(mmat); /* rows of transpose of matrix */
 47:   nm  = mxGetM(mmat);
 48:   nnz = (mxGetJc(mmat))[nn];
 49:   ii  = mxGetJc(mmat);
 50:   jj  = mxGetIr(mmat);
 51:   n   = (PetscInt)nn;
 52:   m   = (PetscInt)nm;
 53:   nz  = (PetscInt)nnz;

 55:   if (mat->rmap->n < 0 && mat->cmap->n < 0) {
 56:     /* matrix has not yet had its size set */
 57:     PetscCall(MatSetSizes(mat, n, m, PETSC_DETERMINE, PETSC_DETERMINE));
 58:     PetscCall(MatSetUp(mat));
 59:   } else {
 60:     PetscCheck(mat->rmap->n == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change size of PETSc matrix %" PetscInt_FMT " to %" PetscInt_FMT, mat->rmap->n, n);
 61:     PetscCheck(mat->cmap->n == m, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change size of PETSc matrix %" PetscInt_FMT " to %" PetscInt_FMT, mat->cmap->n, m);
 62:   }
 63:   if (nz != aij->nz) {
 64:     /* number of nonzeros in matrix has changed, so need new data structure */
 65:     PetscCall(MatSeqXAIJFreeAIJ(mat, &aij->a, &aij->j, &aij->i));
 66:     aij->nz = nz;
 67:     PetscCall(PetscMalloc3(aij->nz, &aij->a, aij->nz, &aij->j, mat->rmap->n + 1, &aij->i));

 69:     aij->singlemalloc = PETSC_TRUE;
 70:   }

 72:   PetscCall(PetscArraycpy(aij->a, mxGetPr(mmat), aij->nz));
 73:   /* MATLAB stores by column, not row so we pass in the transpose of the matrix */
 74:   i = aij->i;
 75:   for (k = 0; k < n + 1; k++) i[k] = (PetscInt)ii[k];
 76:   j = aij->j;
 77:   for (k = 0; k < nz; k++) j[k] = (PetscInt)jj[k];

 79:   for (k = 0; k < mat->rmap->n; k++) aij->ilen[k] = aij->imax[k] = aij->i[k + 1] - aij->i[k];

 81:   mat->nonzerostate++; /* since the nonzero structure can change anytime force the Inode information to always be rebuilt */
 82:   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
 83:   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject obj, void *mengine)
 88: {
 89:   Mat      mat = (Mat)obj;
 90:   mxArray *mmat;

 92:   PetscFunctionBegin;
 93:   mmat = engGetVariable((Engine *)mengine, obj->name);
 94:   PetscCall(MatSeqAIJFromMatlab(mmat, mat));
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: static PetscErrorCode MatSolve_Matlab(Mat A, Vec b, Vec x)
 99: {
100:   const char *_A, *_b, *_x;

102:   PetscFunctionBegin;
103:   /* make sure objects have names; use default if not */
104:   PetscCall(PetscObjectName((PetscObject)b));
105:   PetscCall(PetscObjectName((PetscObject)x));

107:   PetscCall(PetscObjectGetName((PetscObject)A, &_A));
108:   PetscCall(PetscObjectGetName((PetscObject)b, &_b));
109:   PetscCall(PetscObjectGetName((PetscObject)x, &_x));
110:   PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), (PetscObject)b));
111:   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "%s = u%s\\(l%s\\(p%s*%s));", _x, _A, _A, _A, _b));
112:   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "%s = 0;", _b));
113:   /* PetscCall(PetscMatlabEnginePrintOutput(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)),stdout));  */
114:   PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), (PetscObject)x));
115:   PetscFunctionReturn(PETSC_SUCCESS);
116: }

118: static PetscErrorCode MatLUFactorNumeric_Matlab(Mat F, Mat A, const MatFactorInfo *info)
119: {
120:   size_t    len;
121:   char     *_A, *name;
122:   PetscReal dtcol = info->dtcol;

124:   PetscFunctionBegin;
125:   if (F->factortype == MAT_FACTOR_ILU || info->dt > 0) {
126:     /* the ILU form is not currently registered */
127:     if (info->dtcol == PETSC_DEFAULT) dtcol = .01;
128:     F->ops->solve = MatSolve_Matlab;
129:     F->factortype = MAT_FACTOR_LU;

131:     PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), (PetscObject)A));
132:     _A = ((PetscObject)A)->name;
133:     PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "info_%s = struct('droptol',%g,'thresh',%g);", _A, info->dt, dtcol));
134:     PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "[l_%s,u_%s,p_%s] = luinc(%s',info_%s);", _A, _A, _A, _A, _A));
135:     PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "%s = 0;", _A));

137:     PetscCall(PetscStrlen(_A, &len));
138:     PetscCall(PetscMalloc1(len + 2, &name));
139:     PetscCall(PetscSNPrintf(name, len + 2, "_%s", _A));
140:     PetscCall(PetscObjectSetName((PetscObject)F, name));
141:     PetscCall(PetscFree(name));
142:   } else {
143:     PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), (PetscObject)A));
144:     _A = ((PetscObject)A)->name;
145:     PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "[l_%s,u_%s,p_%s] = lu(%s',%g);", _A, _A, _A, _A, dtcol));
146:     PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "%s = 0;", _A));
147:     PetscCall(PetscStrlen(_A, &len));
148:     PetscCall(PetscMalloc1(len + 2, &name));
149:     PetscCall(PetscSNPrintf(name, len + 2, "_%s", _A));
150:     PetscCall(PetscObjectSetName((PetscObject)F, name));
151:     PetscCall(PetscFree(name));

153:     F->ops->solve = MatSolve_Matlab;
154:   }
155:   PetscFunctionReturn(PETSC_SUCCESS);
156: }

158: static PetscErrorCode MatLUFactorSymbolic_Matlab(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
159: {
160:   PetscFunctionBegin;
161:   PetscCheck(A->cmap->N == A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matrix must be square");
162:   F->ops->lufactornumeric = MatLUFactorNumeric_Matlab;
163:   F->assembled            = PETSC_TRUE;
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: static PetscErrorCode MatFactorGetSolverType_seqaij_matlab(Mat A, MatSolverType *type)
168: {
169:   PetscFunctionBegin;
170:   *type = MATSOLVERMATLAB;
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: static PetscErrorCode MatDestroy_matlab(Mat A)
175: {
176:   const char *_A;

178:   PetscFunctionBegin;
179:   PetscCall(PetscObjectGetName((PetscObject)A, &_A));
180:   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(PetscObjectComm((PetscObject)A)), "delete %s l_%s u_%s;", _A, _A, _A));
181:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
182:   PetscFunctionReturn(PETSC_SUCCESS);
183: }

185: static PetscErrorCode MatGetFactor_seqaij_matlab(Mat A, MatFactorType ftype, Mat *F)
186: {
187:   PetscFunctionBegin;
188:   PetscCheck(A->cmap->N == A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matrix must be square");
189:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), F));
190:   PetscCall(MatSetSizes(*F, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
191:   PetscCall(PetscStrallocpy("matlab", &((PetscObject)*F)->type_name));
192:   PetscCall(MatSetUp(*F));

194:   (*F)->ops->destroy           = MatDestroy_matlab;
195:   (*F)->ops->getinfo           = MatGetInfo_External;
196:   (*F)->trivialsymbolic        = PETSC_TRUE;
197:   (*F)->ops->lufactorsymbolic  = MatLUFactorSymbolic_Matlab;
198:   (*F)->ops->ilufactorsymbolic = MatLUFactorSymbolic_Matlab;

200:   PetscCall(PetscObjectComposeFunction((PetscObject)*F, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_matlab));

202:   (*F)->factortype = ftype;
203:   PetscCall(PetscFree((*F)->solvertype));
204:   PetscCall(PetscStrallocpy(MATSOLVERMATLAB, &(*F)->solvertype));
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Matlab(void)
209: {
210:   PetscFunctionBegin;
211:   PetscCall(MatSolverTypeRegister(MATSOLVERMATLAB, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_matlab));
212:   PetscFunctionReturn(PETSC_SUCCESS);
213: }

215: /*MC
216:   MATSOLVERMATLAB - "matlab" - Providing direct solver LU for `MATSEQAIJ` matrix via the external package MATLAB.

218:   Use `./configure` with the options `--with-matlab` to install PETSc with this capability

220:   Options Database Key:
221: . -pc_factor_mat_solver_type matlab - selects MATLAB to do the sparse factorization

223:   Level: beginner

225: .seealso: [](ch_matrices), `Mat`, `PCLU`, `PCFactorSetMatSolverType()`, `MatSolverType`
226: M*/