Changeset 13295


Ignore:
Timestamp:
09/07/12 12:32:15 (13 years ago)
Author:
utke
Message:

NEW: EDF fos forward implementaion for solverx

File:
1 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/src/c/modules/Solverx/SolverxSeq.cpp

    r13261 r13295  
    6666
    6767#ifdef _HAVE_ADOLC_
    68 int EDF_for_solverx(int n, IssmPDouble *x, int m, IssmPDouble *y){/*{{{*/
    69     if(m*(m+1)!=n)_error_("Stiffness matrix should be square!");
    70     SolverxSeq(y,x, x+m*m, m);
     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
    7170    return 0;
    72 }/*}}}*/
     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
    73114void SolverxSeq(IssmDouble *X,IssmDouble *A,IssmDouble *B,int n, Parameters* parameters){/*{{{*/
    74115        // pack inputs to conform to the EDF-prescribed interface
Note: See TracChangeset for help on using the changeset viewer.