Eötvös Quantum Utilities  v4.9.146
Providing the Horsepowers in the Quantum Realm
dgetDiagInv.c
Go to the documentation of this file.
1 /*======================================================================
2  Copyright (C) 2018 Peter Rakyta, Ph.D.
3 
4  This program is free software: you can redistribute it and/or modify
5  it under the terms of the GNU General Public License as published by
6  the Free Software Foundation, either version 3 of the License, or
7  (at your option) any later version.
8 
9  This program is distributed in the hope that it will be useful,
10  but WITHOUT ANY WARRANTY; without even the implied warranty of
11  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  GNU General Public License for more details.
13 
14  You should have received a copy of the GNU General Public License
15  along with this program. If not, see http://www.gnu.org/licenses/.
16 
17 ======================================================================*/
18 
29 #include "mex.h"
30 #include <stdio.h>
31 #include <omp.h>
32 
33 #ifdef __cplusplus
34 extern "C"
35 #endif
36 
38 #ifdef MIC
39 __declspec( target(mic) ) void getdiaginv_mp_dgetdiaginv_(int *sizeA, int *nonzerosA, double *a, int *ia, int *ja, double *invA, int *error);
40 #else
41 void getdiaginv_mp_dgetdiaginv_(int *sizeA, int *nonzerosA, double *a, int *ia, int *ja, double *invA, int *error);
42 #endif
43 
44 
45 
52 void mexFunction(int nlhs, mxArray *plhs[],
53  int nrhs, const mxArray *prhs[])
54 {
55 /* Array information: */
56  mwSize nonzerosA; /*number of nonzeros elements in A*/
57  mwSize sizeAm;
58  mwSize sizeAn; /*size of the matrix A (rows, and cols) */
59  int sizeAn_fort, nonzerosA_fort; /*fortran compatible matrix dimensions*/
60 
61  mwSize *ia, *ja;
62  int *ia_fort, *ja_fort; /* fortran compatible indices */
63  int error;
64  double* a, *invA_fort;
65 
66 
67 /*-----------------------------------------------------------------------*/
68 /* Check for proper number of arguments. */
69  if(nrhs != 1) {
70  mexErrMsgIdAndTxt ("MATLAB:timestwo:nInput", "One input required.");
71  }
72  else if (nlhs > 1) {
73  mexErrMsgIdAndTxt ("MATLAB:timestwo:nOutput", "Too many output arguments.");
74  }
75 
76 /* Validate inputs */
77 
78 /* Check that the input is a double.*/
79  if ((mxIsDouble(prhs[0]) == 0)) {
80  mexErrMsgIdAndTxt ("MATLAB:timestwo:NonDouble", "Input must be a double.");
81  }
82 
83 /* Check if inputs are real.*/
84  if ((mxIsComplex(prhs[0]) != 0) ) {
85  mexErrMsgIdAndTxt ("MATLAB:convec:NonReal", "Input must be real.");
86  }
87 
88 /* Get the size of the input array. (MATLAB) */
89  nonzerosA = mxGetNzmax(prhs[0]);
90  sizeAm = mxGetM(prhs[0]);
91  sizeAn = mxGetN(prhs[0]);
92 
93  ia = mxGetJc(prhs[0]);
94  ja = mxGetIr(prhs[0]);
95  a = mxGetPr(prhs[0]);
96 
97  ia_fort = (int*)malloc( (sizeAm+1)*sizeof(int) );
98  ja_fort = (int*)malloc( nonzerosA*sizeof(int) );
99 
100  mwSize idx;
101  for( idx=1; idx<=sizeAm+1; idx++ ) {
102  ia_fort[idx-1] = (int)ia[idx-1] + 1;
103  }
104 
105  for( idx=1; idx<=nonzerosA; idx++ ) {
106  ja_fort[idx-1] = (int)ja[idx-1] + 1;
107  }
108 
109 
110 
111  plhs[0] = mxCreateDoubleMatrix(sizeAn, 1, mxREAL); /* invA */
112 
113  invA_fort = (double *)mxGetData(plhs[0]);
114 
115  sizeAn_fort = (int)sizeAn;
116  nonzerosA_fort = (int)nonzerosA;
117 
118 #ifdef MIC
119 #pragma offload target(mic) \
120  in(sizeAn_fort, nonzerosA_fort) \
121  in(a:length(nonzerosA)) \
122  in(ia_fort:length(sizeAn+1)) \
123  in(ja_fort:length(nonzerosA)) \
124  out(invA_fort:length(sizeAn_fort)) \
125  out(error)
126  {
127 #endif
128 
129 #ifdef __INTEL_COMPILER
130  getdiaginv_mp_dgetdiaginv_(&sizeAn_fort, &nonzerosA_fort, a, ia_fort, ja_fort, invA_fort, &error);
131 #else
132  __getdiaginv_MOD_dgetdiaginv(&sizeAn_fort, &nonzerosA_fort, a, ia_fort, ja_fort, invA_fort, &error);
133 #endif
134 
135 #ifdef MIC
136  }
137 #endif
138 
139  free(ia_fort);
140  free(ja_fort);
141 
142 /* Check the error flag .*/
143  if ((error != 0) ) {
144  char output[100];
145  sprintf(output, "PARDISO solver ended with error flag: %i\n", error);
146  mexErrMsgIdAndTxt ("EQuUs:dgetPartialInv:PardisoError", output);
147  }
148 
149 #ifdef DEBUG
150  mexPrintf("Pardiso calculations finished for real matrix\n");
151 #endif
152 
153 }
154 
155 
void getdiaginv_mp_dgetdiaginv_(int *sizeA, int *nonzerosA, double *a, int *ia, int *ja, double *invA, int *error)
For specification see getdiaginv::dgetdiaginv.
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
Gateway routine to call the dgetdiaginv subroutine from the mkl_EQuUs package.
Definition: dgetDiagInv.c:52