void test_content() { rng_type * rng = rng_alloc(MZRAN , INIT_DEFAULT); matrix_type * PC = matrix_alloc( 3 , 10); matrix_type * PC_obs = matrix_alloc( 3 , 1 ); double_vector_type * singular_values = double_vector_alloc(3 , 1); matrix_random_init( PC , rng ); matrix_random_init( PC_obs , rng ); { pca_plot_data_type * data = pca_plot_data_alloc("KEY" , PC , PC_obs, singular_values); for (int i=0; i < matrix_get_rows( PC ); i++) { const pca_plot_vector_type * vector = pca_plot_data_iget_vector( data , i ); test_assert_double_equal( matrix_iget( PC_obs , i , 0) , pca_plot_vector_get_obs_value( vector ) ); test_assert_double_equal( double_vector_iget( singular_values , i), pca_plot_vector_get_singular_value( vector ) ); for (int j=0; j < matrix_get_columns( PC ); j++) test_assert_double_equal( matrix_iget( PC , i , j ) , pca_plot_vector_iget_sim_value( vector , j )); test_assert_int_equal( matrix_get_columns( PC ) , pca_plot_vector_get_size( vector )); } pca_plot_data_free( data ); } double_vector_free( singular_values ); matrix_free( PC ); matrix_free( PC_obs ); }
static void dgemm_debug(const matrix_type *C , const matrix_type *A , const matrix_type * B , bool transA, bool transB) { printf("\nC = [%d , %d]\n",matrix_get_rows( C ) , matrix_get_columns(C)); printf("A: [%d , %d]", matrix_get_rows( A ) , matrix_get_columns(A)); if (transA) printf("^T"); printf("\nB: [%d , %d]", matrix_get_rows( B ) , matrix_get_columns(B)); if (transB) printf("^T"); printf("\n\n"); printf("[%d ,%d] = ",matrix_get_rows( C ) , matrix_get_columns(C)); if (transA) printf("[%d ,%d] x ",matrix_get_rows( A ) , matrix_get_columns(A)); else printf("[%d ,%d] x ",matrix_get_columns( A ) , matrix_get_rows(A)); if (transB) printf("[%d ,%d]\n",matrix_get_rows( B ) , matrix_get_columns(B)); else printf("[%d ,%d]\n",matrix_get_columns( B ) , matrix_get_rows(B)); }
/* Updates matrix A by subtracting matrix B - elementwise. */ void matrix_inplace_sub(matrix_type * A , const matrix_type * B) { if ((A->rows == B->rows) && (A->columns == B->columns)) { int i,j; for (j = 0; j < A->columns; j++) for (i=0; i < A->rows; i++) A->data[ GET_INDEX(A,i,j) ] -= B->data[ GET_INDEX(B,i,j) ]; } else util_abort("%s: size mismatch A:[%d,%d] B:[%d,%d]\n",__func__ , matrix_get_rows(A), matrix_get_columns(A), matrix_get_rows(B), matrix_get_columns(B)); }
/* Currently only used as 'support' function for the matrix_det function. */ static void matrix_dgetrf__( matrix_type * A, int * ipiv, int * info) { int lda = matrix_get_column_stride( A ); int m = matrix_get_rows( A ); int n = matrix_get_columns( A ); dgetrf_( &m , &n , matrix_get_data( A ) , &lda , ipiv , info); }
int matrix_inv( matrix_type * A ) { matrix_lapack_assert_square( A ); { int dgetrf_info; int info; int n = matrix_get_columns( A ); int * ipiv = util_malloc( n * sizeof * ipiv ); matrix_dgetrf__( A , ipiv , &dgetrf_info ); { int lda = matrix_get_column_stride( A ); double * work = util_malloc( sizeof * work ); int work_size; /* First call: determine optimal worksize: */ work_size = -1; dgetri_( &n , matrix_get_data( A ), &lda , ipiv , work , &work_size , &info); if (info == 0) { work_size = (int) work[0]; work = util_realloc( work , sizeof * work * work_size ); dgetri_( &n , matrix_get_data( A ), &lda , ipiv , work , &work_size , &info); } else util_abort("%s: dgetri_ returned info:%d \n",__func__ , info); free( work ); } free( ipiv ); return info; } }
void sqrt_enkf_initX(void * module_data , matrix_type * X , matrix_type * A , matrix_type * S , matrix_type * R , matrix_type * dObs , matrix_type * E , matrix_type *D ) { sqrt_enkf_data_type * data = sqrt_enkf_data_safe_cast( module_data ); { int ncomp = std_enkf_get_subspace_dimension( data->std_data ); double truncation = std_enkf_get_truncation( data->std_data ); int nrobs = matrix_get_rows( S ); int ens_size = matrix_get_columns( S ); int nrmin = util_int_min( ens_size , nrobs); matrix_type * W = matrix_alloc(nrobs , nrmin); double * eig = util_calloc( nrmin , sizeof * eig ); matrix_subtract_row_mean( S ); /* Shift away the mean */ enkf_linalg_lowrankCinv( S , R , W , eig , truncation , ncomp); enkf_linalg_init_sqrtX( X , S , data->randrot , dObs , W , eig , false); matrix_free( W ); free( eig ); enkf_linalg_checkX( X , false ); } }
void matrix_dorgqr(matrix_type * A , double * tau, int num_reflectors) { /* num_reflectors == length of tau. */ int lda = matrix_get_column_stride( A ); int m = matrix_get_rows( A ); int n = matrix_get_columns( A ); double * work = util_malloc(sizeof * work ); int worksize; int info; /* Determine optimal worksize. */ worksize = -1; dorgqr_(&m , &n , &num_reflectors , matrix_get_data( A ), &lda , tau , work , &worksize , &info); if (info != 0) util_abort("%s: dorgqf routine failed with info:%d \n",__func__ , info); worksize = ( int ) work[0]; { double * tmp = realloc(work , sizeof * work * worksize ); if (tmp == NULL) { /* OK - we could not get the optimal worksize, try again with the minimum. */ worksize = n; work = util_realloc(work , sizeof * work * worksize ); } else work = tmp; /* The request for optimal worksize succeeded */ } /* Second call - do the actual computation. */ dorgqr_(&m , &n , &num_reflectors , matrix_get_data( A ), &lda , tau , work , &worksize , &info); if (info != 0) util_abort("%s: dorqf routine failed with info:%d \n",__func__ , info); free( work ); }
void matrix_column_compressed_memcpy(matrix_type * target, const matrix_type * src, const bool_vector_type * mask) { if (bool_vector_count_equal( mask , true ) != matrix_get_columns( target )) util_abort("%s: size mismatch. \n",__func__); if (bool_vector_size( mask ) != matrix_get_columns( src)) util_abort("%s: size mismatch. \n",__func__); { int target_col = 0; int src_col; for (src_col = 0; src_col < bool_vector_size( mask ); src_col++) { if (bool_vector_iget( mask , src_col)) { matrix_copy_column( target , src , target_col , src_col); target_col++; } } } }
int matrix_gram_schmidt_process(matrix_t *q, matrix_t *p) { int n; matrix_t *tmp; assert(q); assert(p); assert(matrix_get_columns(q) >= matrix_get_columns(p)); assert(matrix_get_rows(q) >= matrix_get_rows(p)); tmp = matrix_new_and_copy(p); n = matrix_self_gram_schmidt_process(tmp, 0, 0, matrix_get_columns(tmp), matrix_get_rows(tmp)); matrix_copy_matrix(q, 0, 0, tmp, 0, 0, n, matrix_get_rows(tmp)); matrix_destroy(tmp); return n; }
void matrix_gram_set( const matrix_type * X , matrix_type * G, bool col) { int G_rows = matrix_get_rows( G ); int G_cols = matrix_get_columns( G ); int X_rows = matrix_get_rows( X ); int X_cols = matrix_get_columns( X ); if (col) { // Calculate X' · X if ((G_rows == G_cols) && (X_cols == G_rows)) matrix_dgemm( G , X , X , true , false , 1 , 0); else util_abort("%s: dimension mismatch \n",__func__); } else { // Calculate X · X' if ((G_rows == G_cols) && (X_rows == G_rows)) matrix_dgemm( G , X , X , false , true , 1 , 0); else util_abort("%s: dimension mismatch \n",__func__); } }
matrix_t *matrix_new_and_gram_schmidt_process(matrix_t *p) { int n; matrix_t *tmp, *q; assert(p); tmp = matrix_new_and_copy(p); n = matrix_self_gram_schmidt_process(tmp, 0, 0, matrix_get_columns(tmp), matrix_get_rows(tmp)); if (n == matrix_get_columns(tmp)) return tmp; q = matrix_new(n, matrix_get_rows(tmp), false); matrix_copy_matrix(q, 0, 0, tmp, 0, 0, n, matrix_get_rows(tmp)); matrix_destroy(tmp); return q; }
matrix_type * matrix_alloc_column_compressed_copy(const matrix_type * src, const bool_vector_type * mask) { if (bool_vector_size( mask ) != matrix_get_columns( src )) util_abort("%s: size mismatch. Src matrix has %d rows mask has:%d elements\n", __func__ , matrix_get_rows( src ) , bool_vector_size( mask )); { int target_columns = bool_vector_count_equal( mask , true ); matrix_type * target = matrix_alloc( matrix_get_rows( src ) , target_columns ); matrix_column_compressed_memcpy( target , src , mask ); return target; } }
matrix_t *cmatrix_new_and_gram_schmidt_process(matrix_t *p) { int n; matrix_t *q, *tmp; assert(p); assert(matrix_is_imaginary(p)); tmp = matrix_new_and_copy(p); n = cmatrix_self_gram_schmidt_process(tmp, 0, 0, matrix_get_columns(tmp), matrix_get_rows(tmp)); if (n == matrix_get_columns(tmp)) return tmp; q = matrix_new(n, matrix_get_rows(tmp), true); cmatrix_copy_cmatrix(q, 0, 0, tmp, 0, 0, n, matrix_get_rows(tmp)); matrix_destroy(tmp); return q; }
void test_readwrite() { test_work_area_type * test_area = test_work_area_alloc("matrix-test"); { rng_type * rng = rng_alloc(MZRAN , INIT_DEV_URANDOM ); matrix_type * m1 = matrix_alloc(3 , 3); matrix_type * m2 = matrix_alloc(3 , 3); matrix_random_init( m1 , rng ); matrix_assign(m2 , m1); test_assert_true( matrix_equal( m1 , m2 ) ); { FILE * stream = util_fopen("m1" , "w"); matrix_fwrite( m1 , stream ); fclose( stream ); } matrix_random_init( m1 , rng ); test_assert_false( matrix_equal( m1 , m2 ) ); { FILE * stream = util_fopen("m1" , "r"); matrix_free( m1 ); m1 = matrix_alloc(1,1); printf("-----------------------------------------------------------------\n"); matrix_fread( m1 , stream ); test_assert_int_equal( matrix_get_rows(m1) , matrix_get_rows( m2)); test_assert_int_equal( matrix_get_columns(m1) , matrix_get_columns( m2)); util_fseek( stream , 0 , SEEK_SET); { matrix_type * m3 = matrix_fread_alloc( stream ); test_assert_true( matrix_equal( m2 , m3 )); matrix_free( m3 ); } fclose( stream ); } test_assert_true( matrix_equal( m1 , m2 ) ); matrix_free( m2 ); matrix_free( m1 ); rng_free( rng ); } test_work_area_free( test_area ); }
static void * matrix_inplace_matmul_mt__(void * arg) { arg_pack_type * arg_pack = arg_pack_safe_cast( arg ); int row_offset = arg_pack_iget_int( arg_pack , 0 ); int rows = arg_pack_iget_int( arg_pack , 1 ); matrix_type * A = arg_pack_iget_ptr( arg_pack , 2 ); const matrix_type * B = arg_pack_iget_const_ptr( arg_pack , 3 ); matrix_type * A_view = matrix_alloc_shared( A , row_offset , 0 , rows , matrix_get_columns( A )); matrix_inplace_matmul( A_view , B ); matrix_free( A_view ); return NULL; }
void sqrt_enkf_init_update( void * arg , const matrix_type * S , const matrix_type * R , const matrix_type * dObs , const matrix_type * E , const matrix_type * D ) { sqrt_enkf_data_type * sqrt_data = sqrt_enkf_data_safe_cast( arg ); { int ens_size = matrix_get_columns( S ); sqrt_data->randrot = enkf_linalg_alloc_mp_randrot( ens_size , sqrt_data->rng ); } }
matrix_type * matrix_alloc_gram( const matrix_type * X , bool col) { int X_rows = matrix_get_rows( X ); int X_columns = matrix_get_columns( X ); matrix_type * G; if (col) G = matrix_alloc( X_columns , X_columns ); else G = matrix_alloc( X_rows , X_rows ); matrix_gram_set( X , G , col); return G; }
void matrix_dgemv(const matrix_type * A , const double *x , double * y, bool transA , double alpha , double beta) { int m = matrix_get_rows( A ); int n = matrix_get_columns( A ); int lda = matrix_get_column_stride( A ); int incx = 1; int incy = 1; char transA_c; if (transA) transA_c = 'T'; else transA_c = 'N'; dgemv_(&transA_c , &m , &n , &alpha , matrix_get_data( A ) , &lda , x , &incx , &beta , y , &incy); }
void matrix_dgesv(matrix_type * A , matrix_type * B) { matrix_lapack_assert_square( A ); matrix_lapack_assert_fortran_layout( B ); { int n = matrix_get_rows( A ); int lda = matrix_get_column_stride( A ); int ldb = matrix_get_column_stride( B ); int nrhs = matrix_get_columns( B ); long int * ipivot = util_calloc( n , sizeof * ipivot ); int info; dgesv_(&n , &nrhs , matrix_get_data( A ) , &lda , ipivot , matrix_get_data( B ), &ldb , &info); if (info != 0) util_abort("%s: low level lapack routine: dgesv() failed with info:%d \n",__func__ , info); free(ipivot); } }
static void lars_estimate_init( lars_type * lars, matrix_type * X , matrix_type * Y) { int nvar = matrix_get_columns( lars->X ); matrix_assign( X , lars->X ); matrix_assign( Y , lars->Y ); if (lars->X_norm != NULL) matrix_free( lars->X_norm ); lars->X_norm = matrix_alloc(1 , nvar ); if (lars->X_mean != NULL) matrix_free( lars->X_mean ); lars->X_mean = matrix_alloc(1 , nvar ); if (lars->beta != NULL) matrix_free( lars->beta ); lars->beta = matrix_alloc( nvar , nvar ); lars->Y_mean = regression_scale( X , Y , lars->X_mean , lars->X_norm); }
void test_state() { rng_type * rng = rng_alloc( MZRAN , INIT_DEFAULT ); int ens_size = 10; int active_size = 8; int rows = 100; matrix_type * state = matrix_alloc(1,1); bool_vector_type * ens_mask = bool_vector_alloc(ens_size , false); matrix_type * A = matrix_alloc( rows , active_size); matrix_type * A2 = matrix_alloc( rows, active_size ); matrix_type * A3 = matrix_alloc( 1,1 ); for (int i=0; i < active_size; i++) bool_vector_iset( ens_mask , i + 1 , true ); matrix_random_init(A , rng); rml_enkf_common_store_state( state , A , ens_mask ); test_assert_int_equal( matrix_get_rows( state ) , rows ); test_assert_int_equal( matrix_get_columns( state ) , ens_size ); { int g; int a = 0; for (g=0; g < ens_size; g++) { if (bool_vector_iget( ens_mask , g )) { test_assert_true( matrix_columns_equal( state , g , A , a )); a++; } } } rml_enkf_common_recover_state( state , A2 , ens_mask); rml_enkf_common_recover_state( state , A3 , ens_mask); test_assert_true( matrix_equal( A , A2 )); test_assert_true( matrix_equal( A , A3 )); bool_vector_free( ens_mask ); matrix_free( state ); matrix_free( A ); }
double matrix_det( matrix_type *A ) { matrix_lapack_assert_square( A ); { int dgetrf_info; double det = 1; double det_scale = 0; int n = matrix_get_columns( A ); int * ipiv = util_malloc( n * sizeof * ipiv ); matrix_dgetrf__( A , ipiv , &dgetrf_info ); { int i; for (i=0; i < n; i++) { det *= matrix_iget(A , i , i); if (det == 0) return 0; /* Holy f**k - a float == comparison ?? */ if (ipiv[i] != (i + 1)) /* A permutation has taken place. */ det *= -1; /* Try to avoid overflow/underflow by factoring out the order of magnitude. */ while (fabs(det) > 10.0) { det /= 10; det_scale += 1; } while (fabs(det) < 1.0) { det *= 10; det_scale -= 1; } } } free( ipiv ); return det * pow(10 , det_scale ); } }
void matrix_dgemm(matrix_type *C , const matrix_type *A , const matrix_type * B , bool transA, bool transB , double alpha , double beta) { int m = matrix_get_rows( C ); int n = matrix_get_columns( C ); int lda = matrix_get_column_stride( A ); int ldb = matrix_get_column_stride( B ); int ldc = matrix_get_column_stride( C ); char transA_c; char transB_c; int k , innerA, innerB , outerA , outerB; if (transA) k = matrix_get_rows( A ); else k = matrix_get_columns( A ); if (transA) { innerA = matrix_get_rows(A); outerA = matrix_get_columns(A); transA_c = 'T'; } else { innerA = matrix_get_columns(A); outerA = matrix_get_rows(A); transA_c = 'N'; } if (transB) { innerB = matrix_get_columns( B ); outerB = matrix_get_rows( B ); transB_c = 'T'; } else { transB_c = 'N'; innerB = matrix_get_rows( B ); outerB = matrix_get_columns( B ); } /* This is the dimension check which must pass: -------------------------------------------------- A | B | Columns(A) = Rows(B) Trans(A) | Trans(B) | Rows(A) = Columns(B) A | Trans(B) | Columns(A) = Columns(B) Trans(A) | B | Rows(A) = Rows(B) -------------------------------------------------- -------------------------------------------------- A | Rows(A) = Rows(C) Trans(A) | Columns(A) = Rows(C) B | Columns(B) = Columns(C) Trans(B) | Rows(B) = Columns(B) -------------------------------------------------- */ if (innerA != innerB) { dgemm_debug(C,A,B,transA , transB); util_abort("%s: matrix size mismatch between A and B \n", __func__); } if (outerA != matrix_get_rows( C )) { dgemm_debug(C,A,B,transA , transB); printf("outerA:%d rows(C):%d \n",outerA , matrix_get_rows( C )); util_abort("%s: matrix size mismatch between A and C \n",__func__); } if (outerB != matrix_get_columns( C )) { dgemm_debug(C,A,B,transA , transB); util_abort("%s: matrix size mismatch between B and C \n",__func__); } if (!ldc >= util_int_max(1 , m)) { dgemm_debug(C,A,B,transA , transB); fprintf(stderr,"Tried to capture blas message: \"** On entry to DGEMM parameter 13 had an illegal value\"\n"); fprintf(stderr,"m:%d ldc:%d ldc should be >= max(1,%d) \n",m,ldc,m); util_abort("%s: invalid value for ldc\n",__func__); } dgemm_(&transA_c , // 1 &transB_c , // 2 &m , // 3 &n , // 4 &k , // 5 &alpha , // 6 matrix_get_data( A ) , // 7 &lda , // 8 matrix_get_data( B ) , // 9 &ldb , // 10 &beta , // 11 matrix_get_data( C ) , // 12 &ldc); // 13 }
matrix_type * matrix_alloc_transpose( const matrix_type * A) { matrix_type * B = matrix_alloc( matrix_get_columns( A ) , matrix_get_rows( A )); matrix_transpose( A , B ); return B; }
void matrix_inplace_matmul(matrix_type * A, const matrix_type * B) { if ((A->columns == B->rows) && (B->rows == B->columns)) { double * tmp = util_malloc( sizeof * A->data * A->columns ); int i,j,k; for (i=0; i < A->rows; i++) { /* Clearing the tmp vector */ for (k=0; k < B->rows; k++) tmp[k] = 0; for (j=0; j < B->rows; j++) { double scalar_product = 0; for (k=0; k < A->columns; k++) scalar_product += A->data[ GET_INDEX(A,i,k) ] * B->data[ GET_INDEX(B,k,j) ]; /* Assign first to tmp[j] */ tmp[j] = scalar_product; } for (j=0; j < A->columns; j++) A->data[ GET_INDEX(A , i, j) ] = tmp[j]; } free(tmp); } else util_abort("%s: size mismatch: A:[%d,%d] B:[%d,%d]\n",__func__ , matrix_get_rows(A) , matrix_get_columns(A) , matrix_get_rows(B) , matrix_get_columns(B)); }
void bootstrap_enkf_updateA(void * module_data , matrix_type * A , matrix_type * S , matrix_type * R , matrix_type * dObs , matrix_type * E , matrix_type * D ) { bootstrap_enkf_data_type * bootstrap_data = bootstrap_enkf_data_safe_cast( module_data ); { const int num_cpu_threads = 4; int ens_size = matrix_get_columns( A ); matrix_type * X = matrix_alloc( ens_size , ens_size ); matrix_type * A0 = matrix_alloc_copy( A ); matrix_type * S_resampled = matrix_alloc_copy( S ); matrix_type * A_resampled = matrix_alloc( matrix_get_rows(A0) , matrix_get_columns( A0 )); int ** iens_resample = alloc_iens_resample( bootstrap_data->rng , ens_size ); { int ensemble_members_loop; for ( ensemble_members_loop = 0; ensemble_members_loop < ens_size; ensemble_members_loop++) { int unique_bootstrap_components; int ensemble_counter; /* Resample A and meas_data. Here we are careful to resample the working copy.*/ { { int_vector_type * bootstrap_components = int_vector_alloc( ens_size , 0); for (ensemble_counter = 0; ensemble_counter < ens_size; ensemble_counter++) { int random_column = iens_resample[ ensemble_members_loop][ensemble_counter]; int_vector_iset( bootstrap_components , ensemble_counter , random_column ); matrix_copy_column( A_resampled , A0 , ensemble_counter , random_column ); matrix_copy_column( S_resampled , S , ensemble_counter , random_column ); } int_vector_select_unique( bootstrap_components ); unique_bootstrap_components = int_vector_size( bootstrap_components ); int_vector_free( bootstrap_components ); } if (bootstrap_data->doCV) { const bool_vector_type * ens_mask = NULL; cv_enkf_init_update( bootstrap_data->cv_enkf_data , ens_mask , S_resampled , R , dObs , E , D); cv_enkf_initX( bootstrap_data->cv_enkf_data , X , A_resampled , S_resampled , R , dObs , E , D); } else std_enkf_initX(bootstrap_data->std_enkf_data , X , NULL , S_resampled,R, dObs, E,D ); matrix_inplace_matmul_mt1( A_resampled , X , num_cpu_threads ); matrix_inplace_add( A_resampled , A0 ); matrix_copy_column( A , A_resampled, ensemble_members_loop, ensemble_members_loop); } } } free_iens_resample( iens_resample , ens_size); matrix_free( X ); matrix_free( S_resampled ); matrix_free( A_resampled ); matrix_free( A0 ); } }
/** Will not respect strides - that is considered low level data layout. */ static matrix_type * matrix_alloc_copy__( const matrix_type * src , bool safe_mode) { matrix_type * copy = matrix_alloc__( matrix_get_rows( src ), matrix_get_columns( src ) , safe_mode); if (copy != NULL) matrix_assign(copy , src); return copy; }
double matrix_row_column_dot_product(const matrix_type * m1 , int row1 , const matrix_type * m2 , int col2) { if (m1->columns != m2->rows) util_abort("%s: size mismatch: m1:[%d,%d] m2:[%d,%d] \n",__func__ , matrix_get_rows( m1 ) , matrix_get_columns( m1 ) , matrix_get_rows( m2 ) , matrix_get_columns( m2 )); { int k; double sum = 0; for( k = 0; k < m1->columns; k++) sum += m1->data[ GET_INDEX(m1 , row1 , k) ] * m2->data[ GET_INDEX(m2, k , col2) ]; return sum; } }
void rml_enkf_common_initA__( matrix_type * A , matrix_type * S , matrix_type * Cd , matrix_type * E , matrix_type * D , double truncation, double lamda, matrix_type * Udr, double * Wdr, matrix_type * VdTr) { int nrobs = matrix_get_rows( S ); int ens_size = matrix_get_columns( S ); double a = lamda + 1; matrix_type *tmp = matrix_alloc (nrobs, ens_size); double nsc = 1/sqrt(ens_size-1); printf("The lamda Value is %5.5f\n",lamda); printf("The Value of Truncation is %4.2f \n",truncation); matrix_subtract_row_mean( S ); /* Shift away the mean in the ensemble predictions*/ matrix_inplace_diag_sqrt(Cd); matrix_dgemm(tmp, Cd, S,false, false, 1.0, 0.0); matrix_scale(tmp, nsc); printf("The Scaling of data matrix completed !\n "); // SVD(S) = Ud * Wd * Vd(T) int nsign = enkf_linalg_svd_truncation(tmp , truncation , -1 , DGESVD_MIN_RETURN , Wdr , Udr , VdTr); /* After this we only work with the reduced dimension matrices */ printf("The number of siginificant ensembles are %d \n ",nsign); matrix_type * X1 = matrix_alloc( nsign, ens_size); matrix_type * X2 = matrix_alloc (nsign, ens_size ); matrix_type * X3 = matrix_alloc (ens_size, ens_size ); // Compute the matrices X1,X2,X3 and dA enkf_linalg_rml_enkfX1(X1, Udr ,D ,Cd ); //X1 = Ud(T)*Cd(-1/2)*D -- D= -(dk-d0) enkf_linalg_rml_enkfX2(X2, Wdr ,X1 ,a, nsign); //X2 = ((a*Ipd)+Wd^2)^-1 * X1 matrix_free(X1); enkf_linalg_rml_enkfX3(X3, VdTr ,Wdr,X2, nsign); //X3 = Vd *Wd*X2 printf("The X3 matrix is computed !\n "); matrix_type *dA1= matrix_alloc (matrix_get_rows(A), ens_size); matrix_type * Dm = matrix_alloc_copy( A ); matrix_subtract_row_mean( Dm ); /* Remove the mean from the ensemble of model parameters*/ matrix_scale(Dm, nsc); enkf_linalg_rml_enkfdA(dA1, Dm, X3); //dA = Dm * X3 matrix_inplace_add(A,dA1); //dA matrix_free(X3); matrix_free(Dm); matrix_free(dA1); }
matrix_type * matrix_alloc_matmul(const matrix_type * A, const matrix_type * B) { matrix_type * C = matrix_alloc( matrix_get_rows( A ) , matrix_get_columns( B )); matrix_matmul( C , A , B ); return C; }