/*!\file SolverxSeq * \brief implementation of sequential solver using the GSL librarie */ #ifdef HAVE_CONFIG_H #include #else #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!" #endif #include #include "./Solverx.h" #include "../../shared/shared.h" #include "../../include/include.h" #include "../../io/io.h" #ifdef _HAVE_GSL_ #include #endif #ifdef _HAVE_ADOLC_ #include "../../shared/Numerics/adolc_edf.h" #endif void SolverxSeq(SeqVec** puf,SeqMat* Kff, SeqVec* pf, Parameters* parameters){/*{{{*/ #ifdef _HAVE_GSL_ /*Intermediary: */ int M,N,N2,s; SeqVec *uf = NULL; Kff->GetSize(&M,&N); pf->GetSize(&N2); if(N!=N2)_error_("Right hand side vector of size " << N2 << ", when matrix is of size " << M << "-" << N << " !"); if(M!=N)_error_("Stiffness matrix should be square!"); IssmDouble *x = xNew(N); #ifdef _HAVE_ADOLC_ SolverxSeq(x,Kff->matrix,pf->vector,N,parameters); #else SolverxSeq(x,Kff->matrix,pf->vector,N); #endif uf=new SeqVec(x,N); xDelete(x); /*Assign output pointers:*/ *puf=uf; #else _error_("GSL support not compiled in!"); #endif }/*}}}*/ void SolverxSeq(IssmPDouble **pX, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/ /*Allocate output*/ double* X = xNew(n); /*Solve*/ SolverxSeq(X,A,B,n); /*Assign output pointer*/ *pX=X; } /*}}}*/ #ifdef _HAVE_ADOLC_ int EDF_for_solverx(int n, IssmPDouble *x, int m, IssmPDouble *y){ /*{{{*/ SolverxSeq(y,x, x+m*m, m); // x is where the matrix starts, x+m*m is where the right-hand side starts return 0; } /*}}}*/ int EDF_fos_forward_for_solverx(int n, IssmPDouble *inVal, IssmPDouble *inDeriv, int m, IssmPDouble *outVal, IssmPDouble *outDeriv) { /*{{{*/ #ifdef _HAVE_GSL_ // the matrix will be modified by LU decomposition. Use gsl_A copy for (int i=0; i