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