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

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

CHG add in temporary debug output

File size: 7.2 KB
Line 
1/*!\file SolverxSeq
2 * \brief implementation of sequential solver using the GSL librarie
3 */
4
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
10#include <cstring>
11
12#include "./Solverx.h"
13#include "../../shared/shared.h"
14#include "../../include/include.h"
15#include "../../io/io.h"
16
17#ifdef _HAVE_GSL_
18#include <gsl/gsl_linalg.h>
19#endif
20
21#ifdef _HAVE_ADOLC_
22#include "../../shared/Numerics/adolc_edf.h"
23#endif
24
25void SolverxSeq(SeqVec<IssmDouble>** puf,SeqMat<IssmDouble>* Kff, SeqVec<IssmDouble>* pf, Parameters* parameters){/*{{{*/
26
27 #ifdef _HAVE_GSL_
28 /*Intermediary: */
29 int M,N,N2,s;
30 SeqVec<IssmDouble> *uf = NULL;
31
32 Kff->GetSize(&M,&N);
33 pf->GetSize(&N2);
34
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!");
37 IssmDouble *x = xNew<IssmDouble>(N);
38#ifdef _HAVE_ADOLC_
39 SolverxSeq(x,Kff->matrix,pf->vector,N,parameters);
40#else
41 SolverxSeq(x,Kff->matrix,pf->vector,N);
42#endif
43 uf=new SeqVec<IssmDouble>(x,N);
44 xDelete(x);
45
46 /*Assign output pointers:*/
47 *puf=uf;
48
49 #else
50 _error_("GSL support not compiled in!");
51 #endif
52
53}/*}}}*/
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
67#ifdef _HAVE_ADOLC_
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
70 return 0;
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 for (int i=0; i<m*m; ++i)
78 std::cout << "EDF_fos_forward_for_solverx A["<< i << "]=" << inVal[i] << std::endl;
79 for (int i=0; i<m; ++i)
80 std::cout << "EDF_fos_forward_for_solverx b["<< i << "]=" << inVal[i+m*m] << std::endl;
81 double* Acopy = xNew<double>(m*m);
82 xMemCpy(Acopy,inVal,m*m);
83 /*Initialize gsl matrices and vectors: */
84 gsl_matrix_view gsl_A = gsl_matrix_view_array (Acopy,m,m);
85 gsl_vector_view gsl_b = gsl_vector_view_array (inVal+m*m,m); // the right hand side starts at address inVal+m*m
86 gsl_permutation *perm_p = gsl_permutation_alloc (m);
87 int signPerm;
88 // factorize
89 gsl_linalg_LU_decomp (&gsl_A.matrix, perm_p, &signPerm);
90 gsl_vector *gsl_x_p = gsl_vector_alloc (m);
91 // solve for the value
92 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_b.vector, gsl_x_p);
93 /*Copy result*/
94 xMemCpy(outVal,gsl_vector_ptr(gsl_x_p,0),m);
95 gsl_vector_free(gsl_x_p);
96 // solve for the derivatives acc. to A * dx = r with r=db - dA * x
97 // compute the RHS
98 double* r=xNew<double>(m);
99 for (int i=0; i<m; i++) {
100 r[i]=inDeriv[m*m+i]; // this is db[i]
101 for (int j=0;j<m; j++) {
102 r[i]-=inDeriv[i*n+j]*outVal[j]; // this is dA[i][j]*x[j]
103 }
104 }
105 gsl_vector_view gsl_r=gsl_vector_view_array(r,m);
106 gsl_vector *gsl_dx_p = gsl_vector_alloc(m);
107 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_r.vector, gsl_dx_p);
108 xMemCpy(outDeriv,gsl_vector_ptr(gsl_dx_p,0),m);
109 gsl_vector_free(gsl_dx_p);
110 xDelete(r);
111 gsl_permutation_free(perm_p);
112 xDelete(Acopy);
113 #endif
114 return 0;
115}
116/*}}}*/
117
118int EDF_fov_forward_for_solverx(int n, IssmPDouble *inVal, int directionCount, IssmPDouble **inDeriv, int m, IssmPDouble *outVal, IssmPDouble **outDeriv) { /*{{{*/
119#ifdef _HAVE_GSL_
120 // the matrix will be modified by LU decomposition. Use gsl_A copy
121 double* Acopy = xNew<double>(m*m);
122 xMemCpy(Acopy,inVal,m*m);
123 /*Initialize gsl matrices and vectors: */
124 gsl_matrix_view gsl_A = gsl_matrix_view_array (Acopy,m,m);
125 gsl_vector_view gsl_b = gsl_vector_view_array (inVal+m*m,m); // the right hand side starts at address inVal+m*m
126 gsl_permutation *perm_p = gsl_permutation_alloc (m);
127 int signPerm;
128 // factorize
129 gsl_linalg_LU_decomp (&gsl_A.matrix, perm_p, &signPerm);
130 gsl_vector *gsl_x_p = gsl_vector_alloc (m);
131 // solve for the value
132 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_b.vector, gsl_x_p);
133 /*Copy result*/
134 xMemCpy(outVal,gsl_vector_ptr(gsl_x_p,0),m);
135 gsl_vector_free(gsl_x_p);
136 // solve for the derivatives acc. to A * dx = r with r=db - dA * x
137 double* r=xNew<double>(m);
138 gsl_vector *gsl_dx_p = gsl_vector_alloc(m);
139 for (int dir=0;dir<directionCount;++dir) {
140 // compute the RHS
141 for (int i=0; i<m; i++) {
142 r[i]=inDeriv[m*m+i][dir]; // this is db[i]
143 for (int j=0;j<m; j++) {
144 r[i]-=inDeriv[i*n+j][dir]*outVal[j]; // this is dA[i][j]*x[j]
145 }
146 }
147 gsl_vector_view gsl_r=gsl_vector_view_array(r,m);
148 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_r.vector, gsl_dx_p);
149 // reuse r
150 xMemCpy(r,gsl_vector_ptr(gsl_dx_p,0),m);
151 for (int i=0; i<m; i++) {
152 outDeriv[i][dir]=r[i];
153 }
154 }
155 gsl_vector_free(gsl_dx_p);
156 xDelete(r);
157 gsl_permutation_free(perm_p);
158 xDelete(Acopy);
159 #endif
160 return 0;
161}
162/*}}}*/
163
164int EDF_fos_reverse_for_solverx(int m, double *dp_U, int n, double *dp_Z) { /*{{{*/
165 return 0;
166}
167/*}}}*/
168
169void SolverxSeq(IssmDouble *X,IssmDouble *A,IssmDouble *B,int n, Parameters* parameters){/*{{{*/
170 // pack inputs to conform to the EDF-prescribed interface
171 IssmDouble* adoubleEDFin=xNew<IssmDouble>(n*(n+1)); // packed inputs, i.e. matrix and right hand side
172 for(int i=0; i<n*n;i++)adoubleEDFin[i] =A[i]; // pack matrix
173 for(int i=0; i<n; i++)adoubleEDFin[i+n*n]=B[i]; // pack the right hand side
174 IssmPDouble* pdoubleEDFin=xNew<IssmPDouble>(n*(n+1)); // provide space to transfer inputs during call_ext_fct
175 IssmPDouble* pdoubleEDFout=xNew<IssmPDouble>(n); // provide space to transfer outputs during call_ext_fct
176 // call the wrapped solver through the registry entry we retrieve from parameters
177 call_ext_fct(dynamic_cast<GenericParam<Adolc_edf> * >(parameters->FindParamObject(AdolcParamEnum))->GetParameterValue().myEDF_for_solverx_p,
178 n*(n+1), pdoubleEDFin, adoubleEDFin,
179 n, pdoubleEDFout,X);
180 xDelete(adoubleEDFin);
181 xDelete(pdoubleEDFin);
182 xDelete(pdoubleEDFout);
183}
184/*}}}*/
185#endif
186void SolverxSeq(IssmPDouble *X, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/
187#ifdef _HAVE_GSL_
188 /*GSL Matrices and vectors: */
189 int s;
190 gsl_matrix_view a;
191 gsl_vector_view b;
192 gsl_vector *x = NULL;
193 gsl_permutation *p = NULL;
194 /*A will be modified by LU decomposition. Use copy*/
195 for (int i=0; i<n*n; ++i)
196 std::cout << "SolverxSeq A["<< i << "]=" << A[i] << std::endl;
197 for (int i=0; i<n; ++i)
198 std::cout << "SolverxSeq b["<< i << "]=" << B[i] << std::endl;
199 double* Acopy = xNew<double>(n*n);
200 xMemCpy(Acopy,A,n*n);
201
202 /*Initialize gsl matrices and vectors: */
203 a = gsl_matrix_view_array (Acopy,n,n);
204 b = gsl_vector_view_array (B,n);
205 x = gsl_vector_alloc (n);
206
207 /*Run LU and solve: */
208 p = gsl_permutation_alloc (n);
209 gsl_linalg_LU_decomp (&a.matrix, p, &s);
210 gsl_linalg_LU_solve (&a.matrix, p, &b.vector, x);
211
212 //printf ("x = \n");
213 //gsl_vector_fprintf (stdout, x, "%g");
214
215 /*Copy result*/
216 xMemCpy(X,gsl_vector_ptr(x,0),n);
217
218 /*Clean up and assign output pointer*/
219 xDelete(Acopy);
220 gsl_permutation_free(p);
221 gsl_vector_free(x);
222#endif
223}
224/*}}}*/
Note: See TracBrowser for help on using the repository browser.