sparse_hessian_times_B_kronecker_C: remove instances of mexErrMsgTxt

time-shift
Houtan Bastani 2010-09-22 12:51:00 +02:00
parent ca24c92ed7
commit e57056ad14
4 changed files with 43 additions and 38 deletions

View File

@ -528,7 +528,9 @@ zx=[zx; zeros(M_.exo_nbr,np);zeros(M_.exo_det_nbr,np)];
zu=[zu; eye(M_.exo_nbr);zeros(M_.exo_det_nbr,M_.exo_nbr)];
[nrzx,nczx] = size(zx);
rhs = -sparse_hessian_times_B_kronecker_C(hessian,zx);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zx);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
rhs = -rhs;
%lhs
n = M_.endo_nbr+sum(kstate(:,2) > M_.maximum_endo_lag+1 & kstate(:,2) < M_.maximum_endo_lag+M_.maximum_endo_lead+1);
@ -585,7 +587,8 @@ hu = dr.ghu(nstatic+1:nstatic+npred,:);
%kk = kk(1:npred,1:npred);
%rhs = -hessian*kron(zx,zu)-f1*dr.ghxx(end-nyf+1:end,kk(:))*kron(hx(1:npred,:),hu(1:npred,:));
rhs = sparse_hessian_times_B_kronecker_C(hessian,zx,zu);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zx,zu);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
nyf1 = sum(kstate(:,2) == M_.maximum_endo_lag+2);
hu1 = [hu;zeros(np-npred,M_.exo_nbr)];
@ -607,8 +610,8 @@ dr.ghxu = A\rhs;
kk = reshape([1:np*np],np,np);
kk = kk(1:npred,1:npred);
rhs = sparse_hessian_times_B_kronecker_C(hessian,zu);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zu);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
[err, B1] = A_times_B_kronecker_C(B*dr.ghxx,hu1);
mexErrCheck('A_times_B_kronecker_C', err);
@ -655,7 +658,8 @@ for i=1:M_.maximum_endo_lead
[junk,k3a,k3] = ...
find(M_.lead_lag_incidenceordered(M_.maximum_endo_lag+j+1,:));
nk3a = length(k3a);
B1 = sparse_hessian_times_B_kronecker_C(hessian(:,kh(k3,k3)),gu(k3a,:));
[err, B1] = sparse_hessian_times_B_kronecker_C(hessian(:,kh(k3,k3)),gu(k3a,:));
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
RHS = RHS + jacobia_(:,k2)*guu(k2a,:)+B1;
end
% LHS

View File

@ -282,7 +282,9 @@ zx=[zx; zeros(M_.exo_nbr,np);zeros(M_.exo_det_nbr,np)];
zu=[zu; eye(M_.exo_nbr);zeros(M_.exo_det_nbr,M_.exo_nbr)];
[nrzx,nczx] = size(zx);
rhs = -sparse_hessian_times_B_kronecker_C(hessian,zx);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zx);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
rhs = -rhs;
%lhs
n = M_.endo_nbr+sum(kstate(:,2) > M_.maximum_endo_lag+1 & kstate(:,2) < M_.maximum_endo_lag+M_.maximum_endo_lead+1);
@ -339,7 +341,8 @@ hu = dr.ghu(nstatic+1:nstatic+npred,:);
%kk = kk(1:npred,1:npred);
%rhs = -hessian*kron(zx,zu)-f1*dr.ghxx(end-nyf+1:end,kk(:))*kron(hx(1:npred,:),hu(1:npred,:));
rhs = sparse_hessian_times_B_kronecker_C(hessian,zx,zu);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zx,zu);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
nyf1 = sum(kstate(:,2) == M_.maximum_endo_lag+2);
hu1 = [hu;zeros(np-npred,M_.exo_nbr)];
@ -361,8 +364,8 @@ dr.ghxu = A\rhs;
kk = reshape([1:np*np],np,np);
kk = kk(1:npred,1:npred);
rhs = sparse_hessian_times_B_kronecker_C(hessian,zu);
[err, rhs] = sparse_hessian_times_B_kronecker_C(hessian,zu);
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
[err, B1] = A_times_B_kronecker_C(B*dr.ghxx,hu1);
mexErrCheck('A_times_B_kronecker_C', err);
@ -409,7 +412,8 @@ for i=1:M_.maximum_endo_lead
[junk,k3a,k3] = ...
find(M_.lead_lag_incidenceordered(M_.maximum_endo_lag+j+1,:));
nk3a = length(k3a);
B1 = sparse_hessian_times_B_kronecker_C(hessian(:,kh(k3,k3)),gu(k3a,:));
[err, B1] = sparse_hessian_times_B_kronecker_C(hessian(:,kh(k3,k3)),gu(k3a,:));
mexErrCheck('sparse_hessian_times_B_kronecker_C', err);
RHS = RHS + jacobia_(:,k2)*guu(k2a,:)+B1;
end
% LHS

View File

@ -1,14 +1,15 @@
function D = sparse_hessian_times_B_kronecker_C(A,B,C)
%function D = sparse_hessian_times_B_kronecker_C(A,B,C)
function [err, D] = sparse_hessian_times_B_kronecker_C(A,B,C)
%function [err, D] = sparse_hessian_times_B_kronecker_C(A,B,C)
% Computes A * kron(B,C) where A is a sparse matrix.
%
% INPUTS
% A [double] mA*nA matrix.
% B [double] mB*nB matrix.
% C [double] mC*nC matrix.
% A [double] mA*nA matrix.
% B [double] mB*nB matrix.
% C [double] mC*nC matrix.
%
% OUTPUTS
% D [double] mA*(nC*nB) or mA*(nB*nB) matrix.
% err [double] scalar: 1 indicates failure, 0 indicates success
% D [double] mA*(nC*nB) or mA*(nB*nB) matrix.
%
% ALGORITHM
% none.
@ -16,7 +17,7 @@ function D = sparse_hessian_times_B_kronecker_C(A,B,C)
% SPECIAL REQUIREMENTS
% none.
% Copyright (C) 1996-2008 Dynare Team
% Copyright (C) 1996-2010 Dynare Team
%
% This file is part of Dynare.
%
@ -33,6 +34,10 @@ function D = sparse_hessian_times_B_kronecker_C(A,B,C)
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <http://www.gnu.org/licenses/>.
if nargout~=2
error('sparse_hessian_times_B_kronecker_C provides exactly 2 output arguments.')
end
switch nargin
case 3
D = A_times_B_kronecker_C(A,B,C);
@ -41,3 +46,4 @@ switch nargin
otherwise
error('Two or Three input arguments required!')
end
err = 0;

View File

@ -1,5 +1,5 @@
/*
* Copyright (C) 2007-2009 Dynare Team
* Copyright (C) 2007-2010 Dynare Team
*
* This file is part of Dynare.
*
@ -141,18 +141,12 @@ void
mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
{
// Check input and output:
if ((nrhs > 3) || (nrhs < 2))
{
mexErrMsgTxt("Two or Three input arguments required.");
}
if (nlhs > 1)
{
mexErrMsgTxt("Too many output arguments.");
}
if ((nrhs > 3) || (nrhs < 2) || nlhs!=2)
DYN_MEX_FUNC_ERR_MSG_TXT("sparse_hessian_times_B_kronecker_C takes 2 or 3 input arguments and provides exactly 2 output arguments.");
if (!mxIsSparse(prhs[0]))
{
mexErrMsgTxt("First input must be a sparse (dynare) hessian matrix.");
}
DYN_MEX_FUNC_ERR_MSG_TXT("sparse_hessian_times_B_kronecker_C: First input must be a sparse (dynare) hessian matrix.");
// Get & Check dimensions (columns and rows):
mwSize mA, nA, mB, nB, mC, nC;
mA = mxGetM(prhs[0]);
@ -164,16 +158,12 @@ mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
mC = mxGetM(prhs[2]);
nC = mxGetN(prhs[2]);
if (mB*mC != nA)
{
mexErrMsgTxt("Input dimension error!");
}
DYN_MEX_FUNC_ERR_MSG_TXT("Input dimension error!");
}
else // A*kron(B,B) is to be computed.
{
if (mB*mB != nA)
{
mexErrMsgTxt("Input dimension error!");
}
DYN_MEX_FUNC_ERR_MSG_TXT("Input dimension error!");
}
// Get input matrices:
double *B, *C;
@ -190,13 +180,13 @@ mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
double *D;
if (nrhs == 3)
{
plhs[0] = mxCreateDoubleMatrix(mA, nB*nC, mxREAL);
plhs[1] = mxCreateDoubleMatrix(mA, nB*nC, mxREAL);
}
else
{
plhs[0] = mxCreateDoubleMatrix(mA, nB*nB, mxREAL);
plhs[1] = mxCreateDoubleMatrix(mA, nB*nB, mxREAL);
}
D = mxGetPr(plhs[0]);
D = mxGetPr(plhs[1]);
// Computational part:
if (nrhs == 2)
{
@ -206,4 +196,5 @@ mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
{
sparse_hessian_times_B_kronecker_C(isparseA, jsparseA, vsparseA, B, C, D, mA, nA, mB, nB, mC, nC);
}
plhs[0] = mxCreateDoubleScalar(0);
}