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