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

Last change on this file since 13530 was 13530, checked in by Eric.Larour, 12 years ago

CHG and NEW:
macros.h: changed name of ISSMBOOT and ISSMEND macros to IssmBoot and IssmEnd.
AdolcEdf: moved adolc_edf to src/c/classes and renamed it AdolcEdf. Updated all headers and Makefile.am
accordingly.
Noticed dual existence of objects.h in src/c/classes and src/c/classes/objects! Deleted first one and
fixed references accordingly in bamg and kriging.
Created new EnvironmentInit routine to initialize Petsc and MPI. Updated issm.cpp and kriging.cpp accordingly.

File size: 9.3 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 "../../classes/classes.h"
15#include "../../include/include.h"
16#include "../../io/io.h"
17
18#ifdef _HAVE_GSL_
19#include <gsl/gsl_linalg.h>
20#endif
21
22void SolverxSeq(SeqVec<IssmDouble>** puf,SeqMat<IssmDouble>* Kff, SeqVec<IssmDouble>* pf, Parameters* parameters){/*{{{*/
23
24#ifdef _HAVE_GSL_
25 /*Intermediary: */
26 int M,N,N2,s;
27 SeqVec<IssmDouble> *uf = NULL;
28
29 Kff->GetSize(&M,&N);
30 pf->GetSize(&N2);
31
32 if(N!=N2)_error_("Right hand side vector of size " << N2 << ", when matrix is of size " << M << "-" << N << " !");
33 if(M!=N)_error_("Stiffness matrix should be square!");
34#ifdef _HAVE_ADOLC_
35 ensureContiguousLocations(N);
36#endif
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
44 uf=new SeqVec<IssmDouble>(x,N);
45 xDelete(x);
46
47 /*Assign output pointers:*/
48 *puf=uf;
49
50#else
51 _error_("GSL support not compiled in!");
52#endif
53
54}/*}}}*/
55void SolverxSeq(IssmPDouble **pX, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/
56
57 /*Allocate output*/
58 double* X = xNew<double>(n);
59
60 /*Solve*/
61 SolverxSeq(X,A,B,n);
62
63 /*Assign output pointer*/
64 *pX=X;
65}
66/*}}}*/
67
68#ifdef _HAVE_ADOLC_
69int EDF_for_solverx(int n, IssmPDouble *x, int m, IssmPDouble *y){ /*{{{*/
70 SolverxSeq(y,x, x+m*m, m); // x is where the matrix starts, x+m*m is where the right-hand side starts
71 return 0;
72} /*}}}*/
73int EDF_fos_forward_for_solverx(int n, IssmPDouble *inVal, IssmPDouble *inDeriv, int m, IssmPDouble *outVal, IssmPDouble *outDeriv) { /*{{{*/
74#ifdef _HAVE_GSL_
75 // for (int i=0; i<m*m; ++i) std::cout << "EDF_fos_forward_for_solverx A["<< i << "]=" << inVal[i] << std::endl;
76 // for (int i=0; i<m; ++i) std::cout << "EDF_fos_forward_for_solverx b["<< i << "]=" << inVal[i+m*m] << std::endl;
77 // the matrix will be modified by LU decomposition. Use gsl_A copy
78 double* Acopy = xNew<double>(m*m);
79 xMemCpy(Acopy,inVal,m*m);
80 /*Initialize gsl matrices and vectors: */
81 gsl_matrix_view gsl_A = gsl_matrix_view_array (Acopy,m,m);
82 gsl_vector_view gsl_b = gsl_vector_view_array (inVal+m*m,m); // the right hand side starts at address inVal+m*m
83 gsl_permutation *perm_p = gsl_permutation_alloc (m);
84 int signPerm;
85 // factorize
86 gsl_linalg_LU_decomp (&gsl_A.matrix, perm_p, &signPerm);
87 gsl_vector *gsl_x_p = gsl_vector_alloc (m);
88 // solve for the value
89 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_b.vector, gsl_x_p);
90 /*Copy result*/
91 xMemCpy(outVal,gsl_vector_ptr(gsl_x_p,0),m);
92 gsl_vector_free(gsl_x_p);
93 // for (int i=0; i<m; ++i) std::cout << "EDF_fos_forward_for_solverx x["<< i << "]=" << outVal[i] << std::endl;
94 // solve for the derivatives acc. to A * dx = r with r=db - dA * x
95 // compute the RHS
96 double* r=xNew<double>(m);
97 for (int i=0; i<m; i++) {
98 r[i]=inDeriv[m*m+i]; // this is db[i]
99 for (int j=0;j<m; j++) {
100 r[i]-=inDeriv[i*m+j]*outVal[j]; // this is dA[i][j]*x[j]
101 }
102 }
103 gsl_vector_view gsl_r=gsl_vector_view_array(r,m);
104 gsl_vector *gsl_dx_p = gsl_vector_alloc(m);
105 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_r.vector, gsl_dx_p);
106 xMemCpy(outDeriv,gsl_vector_ptr(gsl_dx_p,0),m);
107 gsl_vector_free(gsl_dx_p);
108 xDelete(r);
109 gsl_permutation_free(perm_p);
110 xDelete(Acopy);
111#endif
112 return 0;
113} /*}}}*/
114int EDF_fov_forward_for_solverx(int n, IssmPDouble *inVal, int directionCount, IssmPDouble **inDeriv, int m, IssmPDouble *outVal, IssmPDouble **outDeriv) { /*{{{*/
115#ifdef _HAVE_GSL_
116 // the matrix will be modified by LU decomposition. Use gsl_A copy
117 double* Acopy = xNew<double>(m*m);
118 xMemCpy(Acopy,inVal,m*m);
119 /*Initialize gsl matrices and vectors: */
120 gsl_matrix_view gsl_A = gsl_matrix_view_array (Acopy,m,m);
121 gsl_vector_view gsl_b = gsl_vector_view_array (inVal+m*m,m); // the right hand side starts at address inVal+m*m
122 gsl_permutation *perm_p = gsl_permutation_alloc (m);
123 int signPerm;
124 // factorize
125 gsl_linalg_LU_decomp (&gsl_A.matrix, perm_p, &signPerm);
126 gsl_vector *gsl_x_p = gsl_vector_alloc (m);
127 // solve for the value
128 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_b.vector, gsl_x_p);
129 /*Copy result*/
130 xMemCpy(outVal,gsl_vector_ptr(gsl_x_p,0),m);
131 gsl_vector_free(gsl_x_p);
132 // solve for the derivatives acc. to A * dx = r with r=db - dA * x
133 double* r=xNew<double>(m);
134 gsl_vector *gsl_dx_p = gsl_vector_alloc(m);
135 for (int dir=0;dir<directionCount;++dir) {
136 // compute the RHS
137 for (int i=0; i<m; i++) {
138 r[i]=inDeriv[m*m+i][dir]; // this is db[i]
139 for (int j=0;j<m; j++) {
140 r[i]-=inDeriv[i*m+j][dir]*outVal[j]; // this is dA[i][j]*x[j]
141 }
142 }
143 gsl_vector_view gsl_r=gsl_vector_view_array(r,m);
144 gsl_linalg_LU_solve (&gsl_A.matrix, perm_p, &gsl_r.vector, gsl_dx_p);
145 // reuse r
146 xMemCpy(r,gsl_vector_ptr(gsl_dx_p,0),m);
147 for (int i=0; i<m; i++) {
148 outDeriv[i][dir]=r[i];
149 }
150 }
151 gsl_vector_free(gsl_dx_p);
152 xDelete(r);
153 gsl_permutation_free(perm_p);
154 xDelete(Acopy);
155#endif
156 return 0;
157}
158/*}}}*/
159int EDF_fos_reverse_for_solverx(int m, double *dp_U, int n, double *dp_Z, double* dp_x, double* dp_y) { /*{{{*/
160 // copy to transpose the matrix
161 double* transposed=xNew<double>(m*m);
162 for (int i=0; i<m; ++i) for (int j=0; j<m; ++j) transposed[j*m+i]=dp_x[i*m+j];
163 gsl_matrix_view aTransposed = gsl_matrix_view_array (transposed,m,m);
164 // the adjoint of the solution is our right-hand side
165 gsl_vector_view x_bar=gsl_vector_view_array(dp_U,m);
166 // the last m elements of dp_Z representing the adjoint of the right-hand side we want to compute:
167 gsl_vector_view b_bar=gsl_vector_view_array(dp_Z+m*m,m);
168 gsl_permutation *perm_p = gsl_permutation_alloc (m);
169 int permSign;
170 gsl_linalg_LU_decomp (&aTransposed.matrix, perm_p, &permSign);
171 gsl_linalg_LU_solve (&aTransposed.matrix, perm_p, &x_bar.vector, &b_bar.vector);
172 // now do the adjoint of the matrix
173 for (int i=0; i<m; ++i) for (int j=0; j<m; ++j) dp_Z[i*m+j]-=dp_Z[m*m+i]*dp_y[j];
174 gsl_permutation_free(perm_p);
175 xDelete(transposed);
176 return 0;
177}
178/*}}}*/
179int EDF_fov_reverse_for_solverx(int m, int p, double **dpp_U, int n, double **dpp_Z, double* dp_x, double* dp_y) { /*{{{*/
180 // copy to transpose the matrix
181 double* transposed=xNew<double>(m*m);
182 for (int i=0; i<m; ++i) for (int j=0; j<m; ++j) transposed[j*m+i]=dp_x[i*m+j];
183 gsl_matrix_view aTransposed = gsl_matrix_view_array (transposed,m,m);
184 gsl_permutation *perm_p = gsl_permutation_alloc (m);
185 int permSign;
186 gsl_linalg_LU_decomp (&aTransposed.matrix, perm_p, &permSign);
187 for (int weightsRowIndex=0;weightsRowIndex<p;++weightsRowIndex) {
188 // the adjoint of the solution is our right-hand side
189 gsl_vector_view x_bar=gsl_vector_view_array(dpp_U[weightsRowIndex],m);
190 // the last m elements of dp_Z representing the adjoint of the right-hand side we want to compute:
191 gsl_vector_view b_bar=gsl_vector_view_array(dpp_Z[weightsRowIndex]+m*m,m);
192 gsl_linalg_LU_solve (&aTransposed.matrix, perm_p, &x_bar.vector, &b_bar.vector);
193 // now do the adjoint of the matrix
194 for (int i=0; i<m; ++i) for (int j=0; j<m; ++j) dpp_Z[weightsRowIndex][i*m+j]-=dpp_Z[weightsRowIndex][m*m+i]*dp_y[j];
195 }
196 gsl_permutation_free(perm_p);
197 xDelete(transposed);
198 return 0;
199}
200/*}}}*/
201void SolverxSeq(IssmDouble *X,IssmDouble *A,IssmDouble *B,int n, Parameters* parameters){/*{{{*/
202 // pack inputs to conform to the EDF-prescribed interface
203 // ensure a contiguous block of locations:
204 ensureContiguousLocations(n*(n+1));
205 IssmDouble* adoubleEDFin=xNew<IssmDouble>(n*(n+1)); // packed inputs, i.e. matrix and right hand side
206 for(int i=0; i<n*n;i++)adoubleEDFin[i] =A[i]; // pack matrix
207 for(int i=0; i<n; i++)adoubleEDFin[i+n*n]=B[i]; // pack the right hand side
208 IssmPDouble* pdoubleEDFin=xNew<IssmPDouble>(n*(n+1)); // provide space to transfer inputs during call_ext_fct
209 IssmPDouble* pdoubleEDFout=xNew<IssmPDouble>(n); // provide space to transfer outputs during call_ext_fct
210 // call the wrapped solver through the registry entry we retrieve from parameters
211 call_ext_fct(dynamic_cast<GenericParam<Adolc_edf> * >(parameters->FindParamObject(AdolcParamEnum))->GetParameterValue().myEDF_for_solverx_p,
212 n*(n+1), pdoubleEDFin, adoubleEDFin,
213 n, pdoubleEDFout,X);
214 // for(int i=0; i<n; i++) {ADOLC_DUMP_MACRO(X[i]);}
215 xDelete(adoubleEDFin);
216 xDelete(pdoubleEDFin);
217 xDelete(pdoubleEDFout);
218}
219/*}}}*/
220#endif
221void SolverxSeq(IssmPDouble *X, IssmPDouble *A, IssmPDouble *B,int n){ /*{{{*/
222#ifdef _HAVE_GSL_
223 /*GSL Matrices and vectors: */
224 int s;
225 gsl_matrix_view a;
226 gsl_vector_view b,x;
227 gsl_permutation *p = NULL;
228// for (int i=0; i<n*n; ++i) std::cout << "SolverxSeq A["<< i << "]=" << A[i] << std::endl;
229// for (int i=0; i<n; ++i) std::cout << "SolverxSeq b["<< i << "]=" << B[i] << std::endl;
230 /*A will be modified by LU decomposition. Use copy*/
231 double* Acopy = xNew<double>(n*n);
232 xMemCpy(Acopy,A,n*n);
233
234 /*Initialize gsl matrices and vectors: */
235 a = gsl_matrix_view_array (Acopy,n,n);
236 b = gsl_vector_view_array (B,n);
237 x = gsl_vector_view_array (X,n);
238
239 /*Run LU and solve: */
240 p = gsl_permutation_alloc (n);
241 gsl_linalg_LU_decomp (&a.matrix, p, &s);
242 gsl_linalg_LU_solve (&a.matrix, p, &b.vector, &x.vector);
243
244 /*Clean up and assign output pointer*/
245 xDelete(Acopy);
246 gsl_permutation_free(p);
247#endif
248}
249/*}}}*/
Note: See TracBrowser for help on using the repository browser.