source: issm/trunk-jpl/src/c/modules/Solverx/SolverxSeq.cpp@ 13196

Last change on this file since 13196 was 13196, checked in by utke, 13 years ago

CHG: streamline the allocation of the temporary vector for the solution (single owner allocates & deletes)

File size: 3.3 KB
RevLine 
[12850]1/*!\file SolverxSeq
2 * \brief implementation of sequential solver using the GSL librarie
[11726]3 */
4
[12445]5#ifdef HAVE_CONFIG_H
6 #include <config.h>
7#else
8#error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
9#endif
[12444]10#include <cstring>
11
[11726]12#include "./Solverx.h"
13#include "../../shared/shared.h"
14#include "../../include/include.h"
15#include "../../io/io.h"
16
[12850]17#ifdef _HAVE_GSL_
18#include <gsl/gsl_linalg.h>
19#endif
[11726]20
[13195]21#ifdef _HAVE_ADOLC_
22#include "../../shared/Numerics/adolc_edf.h"
23#endif
[12850]24
[13195]25void SolverxSeq(SeqVec** puf,SeqMat* Kff, SeqVec* pf, Parameters* parameters){/*{{{*/
26
[12850]27 #ifdef _HAVE_GSL_
[11726]28 /*Intermediary: */
29 int M,N,N2,s;
[12450]30 SeqVec *uf = NULL;
[11726]31
32 Kff->GetSize(&M,&N);
33 pf->GetSize(&N2);
34
[13056]35 if(N!=N2)_error_("Right hand side vector of size " << N2 << ", when matrix is of size " << M << "-" << N << " !");
36 if(M!=N)_error_("Stiffness matrix should be square!");
[13196]37 IssmDouble *x = xNew<IssmDouble>(N);
[13195]38#ifdef _HAVE_ADOLC_
[13196]39 SolverxSeq(x,Kff->matrix,pf->vector,N,parameters);
[13195]40#else
[13196]41 SolverxSeq(x,Kff->matrix,pf->vector,N);
[13195]42#endif
[13196]43 uf=new SeqVec(x,N);
44 xDelete(x);
[11726]45
46 /*Assign output pointers:*/
47 *puf=uf;
[12850]48
49 #else
[13056]50 _error_("GSL support not compiled in!");
[12850]51 #endif
52
[12417]53}/*}}}*/
[12988]54#ifdef _HAVE_ADOLC_
[13185]55
56int EDF_for_solverx(int n, IssmPDouble *x, int m, IssmPDouble *y) {
57 if(m*(m+1)!=n)_error_("Stiffness matrix should be square!");
[13196]58 SolverxSeq(y,x, x+m*m, m);
[13185]59 return 0;
60}
61
[13196]62void SolverxSeq(IssmDouble *X,IssmDouble *A,IssmDouble *B,int n, Parameters* parameters){//{{{
[13195]63 // pack inputs to conform to the EDF-prescribed interface
[13196]64 IssmDouble* adoubleEDFin=xNew<IssmDouble>(n*(n+1)); // packed inputs, i.e. matrix and right hand side
65 for(int i=0; i<n*n;i++)adoubleEDFin[i] =A[i]; // pack matrix
66 for(int i=0; i<n; i++)adoubleEDFin[i+n*n]=B[i]; // pack the right hand side
67 IssmPDouble* pdoubleEDFin=xNew<IssmPDouble>(n*(n+1)); // provide space to transfer inputs during call_ext_fct
68 IssmPDouble* pdoubleEDFout=xNew<IssmPDouble>(n); // provide space to transfer outputs during call_ext_fct
[13195]69 // call the wrapped solver through the registry entry we retrieve from parameters
70 call_ext_fct(dynamic_cast<GenericParam<Adolc_edf> * >(parameters->FindParamObject(AdolcParamEnum))->GetParameterValue().myEDF_for_solverx_p,
[13196]71 n*(n+1), pdoubleEDFin, adoubleEDFin,
72 n, pdoubleEDFout,X);
73 xDelete(adoubleEDFin);
74 xDelete(pdoubleEDFin);
75 xDelete(pdoubleEDFout);
[12988]76}
77/*}}}*/
78#endif
[13196]79void SolverxSeq(IssmPDouble * X, IssmPDouble * A, IssmPDouble * B,int n){ //{{{
[12850]80 #ifdef _HAVE_GSL_
[12988]81 /*GSL Matrices and vectors: */
82 int s;
83 gsl_matrix_view a;
84 gsl_vector_view b;
85 gsl_vector *x = NULL;
86 gsl_permutation *p = NULL;
87 /*A will be modified by LU decomposition. Use copy*/
88 double* Acopy = xNew<double>(n*n);
[13196]89 xMemCpy(Acopy,A,n*n);
[12417]90
[12988]91 /*Initialize gsl matrices and vectors: */
92 a = gsl_matrix_view_array (Acopy,n,n);
93 b = gsl_vector_view_array (B,n);
94 x = gsl_vector_alloc (n);
[12417]95
[12988]96 /*Run LU and solve: */
97 p = gsl_permutation_alloc (n);
98 gsl_linalg_LU_decomp (&a.matrix, p, &s);
99 gsl_linalg_LU_solve (&a.matrix, p, &b.vector, x);
[12417]100
[12988]101 //printf ("x = \n");
102 //gsl_vector_fprintf (stdout, x, "%g");
[12417]103
[12988]104 /*Copy result*/
[13196]105 xMemCpy(X,gsl_vector_ptr(x,0),n);
[12417]106
[12988]107 /*Clean up and assign output pointer*/
[13196]108 xDelete(Acopy);
[12988]109 gsl_permutation_free(p);
110 gsl_vector_free(x);
[12850]111 #endif
[12988]112}
113/*}}}*/
Note: See TracBrowser for help on using the repository browser.