booleantype is_square(SUNMatrix A) { if (SUNDenseMatrix_Rows(A) == SUNDenseMatrix_Columns(A)) return SUNTRUE; else return SUNFALSE; }
/*----------------------------------------------------------------- cvDlsDenseDQJac ----------------------------------------------------------------- This routine generates a dense difference quotient approximation to the Jacobian of f(t,y). It assumes that a dense SUNMatrix is stored column-wise, and that elements within each column are contiguous. The address of the jth column of J is obtained via the accessor function SUNDenseMatrix_Column, and this pointer is associated with an N_Vector using the N_VSetArrayPointer function. Finally, the actual computation of the jth column of the Jacobian is done with a call to N_VLinearSum. -----------------------------------------------------------------*/ int cvDlsDenseDQJac(realtype t, N_Vector y, N_Vector fy, SUNMatrix Jac, CVodeMem cv_mem, N_Vector tmp1) { realtype fnorm, minInc, inc, inc_inv, yjsaved, srur; realtype *y_data, *ewt_data; N_Vector ftemp, jthCol; sunindextype j, N; int retval = 0; CVDlsMem cvdls_mem; /* access DlsMem interface structure */ cvdls_mem = (CVDlsMem) cv_mem->cv_lmem; /* access matrix dimension */ N = SUNDenseMatrix_Rows(Jac); /* Rename work vector for readibility */ ftemp = tmp1; /* Create an empty vector for matrix column calculations */ jthCol = N_VCloneEmpty(tmp1); /* Obtain pointers to the data for ewt, y */ ewt_data = N_VGetArrayPointer(cv_mem->cv_ewt); y_data = N_VGetArrayPointer(y); /* Set minimum increment based on uround and norm of f */ srur = SUNRsqrt(cv_mem->cv_uround); fnorm = N_VWrmsNorm(fy, cv_mem->cv_ewt); minInc = (fnorm != ZERO) ? (MIN_INC_MULT * SUNRabs(cv_mem->cv_h) * cv_mem->cv_uround * N * fnorm) : ONE; for (j = 0; j < N; j++) { /* Generate the jth col of J(tn,y) */ N_VSetArrayPointer(SUNDenseMatrix_Column(Jac,j), jthCol); yjsaved = y_data[j]; inc = SUNMAX(srur*SUNRabs(yjsaved), minInc/ewt_data[j]); y_data[j] += inc; retval = cv_mem->cv_f(t, y, ftemp, cv_mem->cv_user_data); cvdls_mem->nfeDQ++; if (retval != 0) break; y_data[j] = yjsaved; inc_inv = ONE/inc; N_VLinearSum(inc_inv, ftemp, -inc_inv, fy, jthCol); /* DENSE_COL(Jac,j) = N_VGetArrayPointer(jthCol); /\*UNNECESSARY?? *\/ */ } /* Destroy jthCol vector */ N_VSetArrayPointer(NULL, jthCol); /* SHOULDN'T BE NEEDED */ N_VDestroy(jthCol); return(retval); }
int SUNLinSolSolve_LapackDense(SUNLinearSolver S, SUNMatrix A, N_Vector x, N_Vector b, realtype tol) { int n, one, ier; realtype *xdata; if ( (A == NULL) || (S == NULL) || (x == NULL) || (b == NULL) ) return(SUNLS_MEM_NULL); /* copy b into x */ N_VScale(ONE, b, x); /* access x data array */ xdata = N_VGetArrayPointer(x); if (xdata == NULL) { LASTFLAG(S) = SUNLS_MEM_FAIL; return(LASTFLAG(S)); } /* Call LAPACK to solve the linear system */ n = SUNDenseMatrix_Rows(A); one = 1; xgetrs_f77("N", &n, &one, SUNDenseMatrix_Data(A), &n, PIVOTS(S), xdata, &n, &ier, 1); LASTFLAG(S) = (long int) ier; if (ier < 0) return(SUNLS_PACKAGE_FAIL_UNREC); LASTFLAG(S) = SUNLS_SUCCESS; return(LASTFLAG(S)); }
CAMLprim value sunml_lsolver_dense(value vnvec, value vdmat) { CAMLparam2(vnvec, vdmat); #if SUNDIALS_LIB_VERSION >= 300 SUNMatrix dmat = MAT_VAL(vdmat); SUNLinearSolver ls = SUNDenseLinearSolver(NVEC_VAL(vnvec), dmat); if (ls == NULL) { if (SUNDenseMatrix_Rows(dmat) != SUNDenseMatrix_Columns(dmat)) caml_raise_constant(LSOLVER_EXN(MatrixNotSquare)); if (SUNDenseMatrix_Rows(dmat) != NV_LENGTH_S(NVEC_VAL(vnvec))) caml_raise_constant(LSOLVER_EXN(MatrixVectorMismatch)); caml_raise_out_of_memory(); } CAMLreturn(alloc_lsolver(ls)); #else CAMLreturn(Val_unit); #endif }
int SUNLinSolSetup_LapackDense(SUNLinearSolver S, SUNMatrix A) { int n, ier; /* check for valid inputs */ if ( (A == NULL) || (S == NULL) ) return(SUNLS_MEM_NULL); /* Ensure that A is a dense matrix */ if (SUNMatGetID(A) != SUNMATRIX_DENSE) { LASTFLAG(S) = SUNLS_ILL_INPUT; return(LASTFLAG(S)); } /* Call LAPACK to do LU factorization of A */ n = SUNDenseMatrix_Rows(A); xgetrf_f77(&n, &n, SUNDenseMatrix_Data(A), &n, PIVOTS(S), &ier); LASTFLAG(S) = (long int) ier; if (ier > 0) return(SUNLS_LUFACT_FAIL); if (ier < 0) return(SUNLS_PACKAGE_FAIL_UNREC); return(SUNLS_SUCCESS); }
SUNLinearSolver SUNLinSol_LapackDense(N_Vector y, SUNMatrix A) { SUNLinearSolver S; SUNLinearSolver_Ops ops; SUNLinearSolverContent_LapackDense content; sunindextype MatrixRows, VecLength; /* Check compatibility with supplied SUNMatrix and N_Vector */ if (SUNMatGetID(A) != SUNMATRIX_DENSE) return(NULL); if (SUNDenseMatrix_Rows(A) != SUNDenseMatrix_Columns(A)) return(NULL); MatrixRows = SUNDenseMatrix_Rows(A); if ( (N_VGetVectorID(y) != SUNDIALS_NVEC_SERIAL) && (N_VGetVectorID(y) != SUNDIALS_NVEC_OPENMP) && (N_VGetVectorID(y) != SUNDIALS_NVEC_PTHREADS) ) return(NULL); /* optimally this function would be replaced with a generic N_Vector routine */ VecLength = GlobalVectorLength_LapDense(y); if (MatrixRows != VecLength) return(NULL); /* Create linear solver */ S = NULL; S = (SUNLinearSolver) malloc(sizeof *S); if (S == NULL) return(NULL); /* Create linear solver operation structure */ ops = NULL; ops = (SUNLinearSolver_Ops) malloc(sizeof(struct _generic_SUNLinearSolver_Ops)); if (ops == NULL) { free(S); return(NULL); } /* Attach operations */ ops->gettype = SUNLinSolGetType_LapackDense; ops->initialize = SUNLinSolInitialize_LapackDense; ops->setup = SUNLinSolSetup_LapackDense; ops->solve = SUNLinSolSolve_LapackDense; ops->lastflag = SUNLinSolLastFlag_LapackDense; ops->space = SUNLinSolSpace_LapackDense; ops->free = SUNLinSolFree_LapackDense; ops->setatimes = NULL; ops->setpreconditioner = NULL; ops->setscalingvectors = NULL; ops->numiters = NULL; ops->resnorm = NULL; ops->resid = NULL; /* Create content */ content = NULL; content = (SUNLinearSolverContent_LapackDense) malloc(sizeof(struct _SUNLinearSolverContent_LapackDense)); if (content == NULL) { free(ops); free(S); return(NULL); } /* Fill content */ content->N = MatrixRows; content->last_flag = 0; content->pivots = NULL; content->pivots = (sunindextype *) malloc(MatrixRows * sizeof(sunindextype)); if (content->pivots == NULL) { free(content); free(ops); free(S); return(NULL); } /* Attach content and ops */ S->content = content; S->ops = ops; return(S); }