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

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

NEW: EDF fos forward implementaion for solverx

File size: 5.1 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
[13216]25void SolverxSeq(SeqVec<IssmDouble>** puf,SeqMat<IssmDouble>* Kff, SeqVec<IssmDouble>* pf, Parameters* parameters){/*{{{*/
[13195]26
[12850]27 #ifdef _HAVE_GSL_
[11726]28 /*Intermediary: */
29 int M,N,N2,s;
[13216]30 SeqVec<IssmDouble> *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
[13216]43 uf=new SeqVec<IssmDouble>(x,N);
[13196]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}/*}}}*/
[13261]54void SolverxSeq(IssmPDouble **pX, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/
55
56 /*Allocate output*/
57 double* X = xNew<double>(n);
58
59 /*Solve*/
60 SolverxSeq(X,A,B,n);
61
62 /*Assign output pointer*/
63 *pX=X;
64}
65/*}}}*/
66
[12988]67#ifdef _HAVE_ADOLC_
[13295]68int EDF_for_solverx(int n, IssmPDouble *x, int m, IssmPDouble *y){ /*{{{*/
69 SolverxSeq(y,x, x+m*m, m); // x is where the matrix starts, x+m*m is where the right-hand side starts
[13185]70 return 0;
[13295]71}
72/*}}}*/
73
74int EDF_fos_forward_for_solverx(int n, IssmPDouble *inVal, IssmPDouble *inDeriv, int m, IssmPDouble *outVal, IssmPDouble *outDeriv) { /*{{{*/
75#ifdef _HAVE_GSL_
76 // the matrix will be modified by LU decomposition. Use gsl_A copy
77 double* Acopy = xNew<double>(m*m);
78 xMemCpy(Acopy,inVal,m*m);
79 /*Initialize gsl matrices and vectors: */
80 gsl_matrix_view gsl_A = gsl_matrix_view_array (Acopy,m,m);
81 gsl_vector_view gsl_b = gsl_vector_view_array (inVal+m*m,m); // the right hand side starts at address inVal+m*m
82 gsl_permutation *perm_p = gsl_permutation_alloc (m);
83 int signPerm;
84 // factorize
85 gsl_linalg_LU_decomp (&gsl_A.matrix, perm_p, &signPerm);
86 gsl_vector *gsl_x_p = gsl_vector_alloc (m);
87 // solve for the value
88 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_b.vector, gsl_x_p);
89 /*Copy result*/
90 xMemCpy(outVal,gsl_vector_ptr(gsl_x_p,0),m);
91 gsl_vector_free(gsl_x_p);
92 // solve for the derivatives acc. to A * dx = r with r=db - dA * x
93 // compute the RHS
94 double* r=xNew<double>(m);
95 for (int i=0; i<m; i++) {
96 r[i]=inDeriv[m*m+i]; // this is db[i]
97 for (int j=0;j<m; j++) {
98 r[i]-=inDeriv[i*n+j]*outVal[j]; // this is dA[i][j]*x[j]
99 }
100 }
101 gsl_vector_view gsl_r=gsl_vector_view_array(r,m);
102 gsl_vector *gsl_dx_p = gsl_vector_alloc(m);
103 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_r.vector, gsl_dx_p);
104 xMemCpy(outDeriv,gsl_vector_ptr(gsl_dx_p,0),m);
105 gsl_vector_free(gsl_dx_p);
106 xDelete(r);
107 gsl_permutation_free(perm_p);
108 xDelete(Acopy);
109 #endif
110 return 0;
111}
112/*}}}*/
113
[13261]114void SolverxSeq(IssmDouble *X,IssmDouble *A,IssmDouble *B,int n, Parameters* parameters){/*{{{*/
[13195]115 // pack inputs to conform to the EDF-prescribed interface
[13261]116 IssmDouble* adoubleEDFin=xNew<IssmDouble>(n*(n+1)); // packed inputs, i.e. matrix and right hand side
117 for(int i=0; i<n*n;i++)adoubleEDFin[i] =A[i]; // pack matrix
118 for(int i=0; i<n; i++)adoubleEDFin[i+n*n]=B[i]; // pack the right hand side
[13196]119 IssmPDouble* pdoubleEDFin=xNew<IssmPDouble>(n*(n+1)); // provide space to transfer inputs during call_ext_fct
[13261]120 IssmPDouble* pdoubleEDFout=xNew<IssmPDouble>(n); // provide space to transfer outputs during call_ext_fct
[13195]121 // call the wrapped solver through the registry entry we retrieve from parameters
122 call_ext_fct(dynamic_cast<GenericParam<Adolc_edf> * >(parameters->FindParamObject(AdolcParamEnum))->GetParameterValue().myEDF_for_solverx_p,
[13196]123 n*(n+1), pdoubleEDFin, adoubleEDFin,
124 n, pdoubleEDFout,X);
125 xDelete(adoubleEDFin);
126 xDelete(pdoubleEDFin);
127 xDelete(pdoubleEDFout);
[12988]128}
129/*}}}*/
130#endif
[13261]131void SolverxSeq(IssmPDouble *X, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/
132#ifdef _HAVE_GSL_
[12988]133 /*GSL Matrices and vectors: */
134 int s;
135 gsl_matrix_view a;
136 gsl_vector_view b;
137 gsl_vector *x = NULL;
138 gsl_permutation *p = NULL;
139 /*A will be modified by LU decomposition. Use copy*/
140 double* Acopy = xNew<double>(n*n);
[13196]141 xMemCpy(Acopy,A,n*n);
[12417]142
[12988]143 /*Initialize gsl matrices and vectors: */
144 a = gsl_matrix_view_array (Acopy,n,n);
145 b = gsl_vector_view_array (B,n);
146 x = gsl_vector_alloc (n);
[12417]147
[12988]148 /*Run LU and solve: */
149 p = gsl_permutation_alloc (n);
150 gsl_linalg_LU_decomp (&a.matrix, p, &s);
151 gsl_linalg_LU_solve (&a.matrix, p, &b.vector, x);
[12417]152
[12988]153 //printf ("x = \n");
154 //gsl_vector_fprintf (stdout, x, "%g");
[12417]155
[12988]156 /*Copy result*/
[13196]157 xMemCpy(X,gsl_vector_ptr(x,0),n);
[12417]158
[12988]159 /*Clean up and assign output pointer*/
[13196]160 xDelete(Acopy);
[12988]161 gsl_permutation_free(p);
162 gsl_vector_free(x);
[13261]163#endif
[12988]164}
165/*}}}*/
Note: See TracBrowser for help on using the repository browser.