Add test functions for nonlinear solvers.

trust-region-mex
Stéphane Adjemian (Ryûk) 2021-07-23 11:56:24 +02:00
parent 0d092a36a0
commit a905539f60
Signed by: stepan
GPG Key ID: 295C1FE89E17EB3C
15 changed files with 851 additions and 1 deletions

View File

@ -1305,7 +1305,21 @@ EXTRA_DIST = \
run_kronecker_tests.m \
kronecker/test_kron.m \
kronecker/nash_matrices.mat \
deterministic_simulations/pfwee.csv
deterministic_simulations/pfwee.csv \
solver-test-functions/brown.m \
solver-test-functions/broydenbanded.m \
solver-test-functions/broydentridiagonal.m \
solver-test-functions/chebyquad.m \
solver-test-functions/discreteboundaryvalue.m \
solver-test-functions/discreteintegralequation.m \
solver-test-functions/helicalvalley.m \
solver-test-functions/powell1.m \
solver-test-functions/powell2.m \
solver-test-functions/rosenbrock.m \
solver-test-functions/trigonometric.m \
solver-test-functions/variablydimensioned.m \
solver-test-functions/watson.m \
solver-test-functions/wood.m
if ENABLE_MATLAB
check-local: check-matlab

View File

@ -0,0 +1,66 @@
function [fval, fjac] = brown(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
sum = -(n+1);
prod = 1;
for j=1:n
sum = sum+x(j);
prod = prod*x(j);
end
for k=1:n-1
fval(k) = x(k)+sum;
end
fval(n) = prod-1;
if nargout>1
fjac = zeros(n);
prod = 1;
for j=1:n
prod = prod*x(j);
for k=1:n
fjac(k,j) = 1;
end
fjac(j,j) = 2;
end
for j=1:n
tmp = x(j);
if abs(tmp)<eps()
tmp = 1;
prod = 1;
for k=1:n
if k~=j
prod = prod*x(k);
end
end
end
fjac(n,j) = prod/tmp;
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = 0.5*ones(length(x), 1);
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,60 @@
function [fval, fjac] = broydenbanded(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
ml = 5;
mu = 1;
for k=1:n
k1 = max(1, k-ml);
k2 = min(k+mu, n);
tmp = 0;
for j=k1:k2
if j~=k
tmp = tmp + x(j)*(1+x(j));
end
end
fval(k) = x(k)*(2+5*x(k)*x(k))+1-tmp;
end
if nargout>1
fjac = zeros(n);
ml = 5;
mu = 1;
for k=1:n
k1 = max(1,k-ml);
k2 = min(k+mu,n);
for j=k1:k2
if j~=k
fjac(k,j) = -(1+2*x(j));
end
end
fjac(k,k) = 2+15*x(k)*x(k);
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = broydentridiagonal(nan(n,1));
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,56 @@
function [fval, fjac] = broydentridiagonal(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
for k=1:n
tmp = (3-2*x(k))*x(k);
tmp1 = 0;
if k~=1
tmp1 = x(k-1);
end
tmp2 = 0;
if k~=n
tmp2 = x(k+1);
end
fval(k) = tmp-tmp1-2*tmp2+1;
end
if nargout>1
fjac = zeros(n);
for k=1:n
fjac(k,k) = 3-4*x(k);
if k~=1
fjac(k,k-1) = -1;
end
if k~=n
fjac(k,k+1) = -2;
end
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = -ones(n, 1);
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,76 @@
function [fval, fjac] = chebyquad(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
for j=1:n
tmp1 = 1;
tmp2 = 2*x(j)-1;
tmp = 2*tmp2;
for i=1:n
fval(i) = fval(i)+tmp2;
ti = tmp*tmp2-tmp1;
tmp1 = tmp2;
tmp2 = ti;
end
end
tk = 1.0/n;
iev = -1;
for k=1:n
fval(k) = fval(k)*tk;
if iev>0
fval(k) = fval(k)+1/(k^2-1);
end
iev = -iev;
end
if nargout>1
fjac = zeros(n);
tk = 1.0/n;
for j=1:n
tmp1 = 1;
tmp2 = 2*x(j)-1;
tmp = 2*tmp2;
tmp3 = 0;
tmp4 = 2;
for k=1:n
fjac(k,j) = tk*tmp4;
ti = 4*tmp2+tmp*tmp4-tmp3;
tmp3 = tmp4;
tmp4 = ti;
ti = tmp*tmp2-tmp1;
tmp1 = tmp2;
tmp2 = ti;
end
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
h = 1.0/(n+1);
fval = zeros(length(x), 1);
for j=1:n
fval(j) = j*h;
end
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,66 @@
function [fval, fjac] = discreteboundaryvalue(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
h = 1.0/(n+1);
hh = h*h;
for k=1:n
tmp = (x(k)+k*h+1)^3;
tmp1 = 0;
if k~=1
tmp1 = x(k-1);
end
tmp2 = 0;
if k~=n
tmp2 = x(k+1);
end
fval(k) = 2*x(k)-tmp1-tmp2+tmp*hh/2;
end
if nargout>1
fjac = zeros(n);
h = 1/(n+1);
hh = h*h;
for k=1:n
tmp = 3*(x(k)+k*h+1)^2;
fjac(k,k) = 2+tmp*hh/2;
if k~=1
fjac(k,k-1) = -1;
end
if k~=n
fjac(k,k+1) = -1;
end
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
h = 1.0/(n+1);
fval = zeros(n, 1);
for j=1:n
tj = j*h;
fval(j) = tj*(tj-1);
end
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,65 @@
function [fval, fjac] = discreteintegralequation(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
h = 1.0/(n+1);
for k=1:n
tk = k*h;
sum1 = 0;
for j=1:k
tj = j*h;
tmp = (x(j)+tj+1)^3;
sum1 = sum1+tj*tmp;
end
sum2 = 0;
kp1 = k+1;
if n>=kp1
for j=kp1:n
tj = j*h;
tmp = (x(j)+tj+1)^3;
sum2 = sum2+(1-tj)*tmp;
end
end
fval(k) = x(k) + h*((1-tk)*sum1+tk*sum2)/2;
end
if nargout>1
fjac = zeros(n);
h = 1.0/(n+1);
for k=1:n
tk = k*h;
for j=1:n
tj = j*h;
tmp = 3*(x(j)+tj+1)^2;
fjac(k,j) = h*min(tj*(1-tk), tk*(1-tj))*tmp/2;
end
fjac(k,k) = fjac(k,k)+1;
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = discreteboundaryvalue(nan(n,1));
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,60 @@
function [fval, fjac] = helicalvalley(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
C7 = 0.25;
C8 = 0.5;
TPI = 8*atan(1);
if nargin==1
fval = zeros(3, 1);
tmp1 = C7*sign(x(2));
if x(1)>0
tmp1 = atan(x(2)/x(1))/TPI;
end
if x(1)<0
tmp1 = atan(x(2)/x(1))/TPI+C8;
end
tmp2 = sqrt(x(1)*x(1)+x(2)*x(2));
fval(1) = 10*(x(3)-10*tmp1);
fval(2) = 10*(tmp2-1);
fval(3) = x(3);
if nargout>1
fjac = zeros(3);
tmp = x(1)*x(1)+x(2)*x(2);
tmp1 = TPI*tmp;
tmp2 = sqrt(tmp);
fjac(1,1) = 100*x(2)/tmp1;
fjac(1,2) = -100*x(1)/tmp1;
fjac(1,3) = 10;
fjac(2,1) = 10*x(1)/tmp2;
fjac(2,2) = 10*x(2)/tmp2;
fjac(2,3) = 0;
fjac(3,1) = 0;
fjac(3,2) = 0;
fjac(3,3) = 1;
end
elseif nargin
if nargout==1
fval = [-1; 0; 0];
else
error('One output is required for initialization mode.')
end
else
error('Wrong number of input arguments.')
end

View File

@ -0,0 +1,47 @@
function [fval, fjac] = powell1(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
if nargin==1
fval = zeros(4, 1);
fval(1) = x(1)+10*x(2);
fval(2) = sqrt(5)*(x(3)-x(4));
tmp = x(2)-2*x(3);
fval(3) = tmp*tmp;
tmp = x(1)-x(4);
fval(4) = sqrt(10)*tmp*tmp;
if nargout>1
fjac = zeros(4);
fjac(1,1) = 1;
fjac(1,2) = 10;
fjac(2,3) = sqrt(5);
fjac(2,4) = -fjac(2,3);
fjac(3,2) = 2*(x(2)-2*x(3));
fjac(3,3) = -2*fjac(3,2);
fjac(4,1) = 2*sqrt(10)*(x(1)-x(4));
fjac(4,4) = -fjac(4,1);
end
elseif nargin
if nargout==1
fval = [3; -1; 0; 1];
else
error('One output is required for initialization mode.')
end
else
error('Wrong number of input arguments.')
end

View File

@ -0,0 +1,39 @@
function [fval, fjac] = powell2(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
if nargin==1
fval = zeros(2, 1);
fval(1) = 10000*x(1)*x(2)-1;
fval(2) = exp(-x(1))+exp(-x(2))-1.0001;
if nargout>1
fjac = zeros(2);
fjac(1,1) = 10000*x(2);
fjac(1,2) = 10000*x(1);
fjac(2,1) = -exp(x(1));
fjac(2,2) = -exp(x(2));
end
elseif nargin
if nargout==1
fval = [0; 1];
else
error('One output is required for initialization mode.')
end
else
error('Wrong number of input arguments.')
end

View File

@ -0,0 +1,39 @@
function [fval, fjac] = rosenbrock(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
if nargin==1
fval = zeros(2, 1);
fval(1) = 1.0-x(1);
fval(2) = 10.0*(x(2)-fval(1)*fval(1));
if nargout>1
fjac = zeros(2);
fjac(1,1) = -1;
fjac(1,2) = 0;
fjac(2,1) = -20*x(1);
fjac(2,2) = 10;
end
elseif nargin
if nargout==1
fval = [10000; 1];
else
error('One output is required for initialization mode.')
end
else
error('Wrong number of input arguments.')
end

View File

@ -0,0 +1,50 @@
function [fval, fjac] = trigonometric(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
sum = 0;
for j=1:n
fval(j) = cos(x(j));
sum = sum+fval(j);
end
for k=1:n
fval(k) = (n+k)-sin(x(k))-sum-k*fval(k);
end
if nargout>1
fjac = zeros(n);
for j=1:n
tmp = sin(x(j));
for k=1:n
fjac(k,j) = tmp;
end
fjac(j,j) = (j+1)*tmp-cos(x(j));
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = (1/n)*ones(n, 1);
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,59 @@
function [fval, fjac] = variablydimensioned(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
sum = 0;
for j=1:n
sum = sum+j*(x(j)-1);
end
tmp = sum*(1+2*sum*sum);
for k=1:n
fval(k) = x(k)-1+k*tmp;
end
if nargout>1
fjac = zeros(n);
sum = 0;
for j=1:n
sum = sum + j*(x(j)-1);
end
tmp = 1+6*sum*sum;
for k=1:n
for j=1:n
fjac(k,j) = k*j*tmp;
fjac(j,k) = fjac(k,j);
end
fjac(k,k) = fjac(k,k)+1;
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
h = 1.0/n;
fval = zeros(n, 1);
for j=1:n
fval(j) = 1-j*h;
end
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,97 @@
function [fval, fjac] = watson(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
n = length(x);
C9 = 29.0;
if isnumeric(x) && isvector(x) && all(~isnan(x))
fval = zeros(n, 1);
for i=1:29
ti = i/C9;
sum1 = 0;
tmp = 1;
for j=2:n
sum1 = sum1+(j-1)*tmp*x(j);
tmp = tmp*ti;
end
sum2 = 0;
tmp = 1;
for j=1:n
sum2 = sum2+tmp*x(j);
tmp = tmp*ti;
end
tmp1 = sum1-sum2*sum2-1;
tmp2 = 2*ti*sum2;
tmp = 1/ti;
for k=1:n
fval(k) = fval(k) + tmp*((k-1)-tmp2)*tmp1;
tmp = tmp*ti;
end
end
tmp = x(2)-x(1)*x(1)-1;
fval(1) = fval(1)+ x(1)*(1-2*tmp);
fval(2) = fval(2)+tmp;
if nargout>1
fjac = zeros(n);
for i=1:29
ti = i/C9;
sum1 = 0;
tmp = 1;
for j=2:n
sum1 = sum1 + (j-1)*tmp*x(j);
tmp = tmp*ti;
end
sum2 = 0;
tmp = 1;
for j=1:n
sum2 = sum2+tmp*x(j);
tmp = tmp*ti;
end
tmp1 = 2*(sum1-sum2*sum2-1);
tmp2 = 2*sum2;
tmp = ti*ti;
tk = 1;
for k=1:n
tj = tk;
for j=k:n
fjac(k,j) = fjac(k,j)+tj*(((k-1)/ti-tmp2)*((j-1)/ti-tmp2)-tmp1);
tj = tj*ti;
end
tk = tk*tmp;
end
end
fjac(1,1) = fjac(1,1)+6*x(1)*x(1)-2*x(2)+3;
fjac(1,2) = fjac(1,2)-2*x(1);
fjac(2,2) = fjac(2,2)+1;
for k=1:n
for j=k:n
fjac(j,k) = fjac(k,j);
end
end
end
elseif isnumeric(x) && isvector(x) && all(isnan(x))
if nargout==1
fval = zeros(length(x), 1);
else
error('One output is required for initialization mode.')
end
else
error('Wrong input argument.')
end

View File

@ -0,0 +1,56 @@
function [fval, fjac] = wood(x)
% Copyright © 2021 Dynare Team
%
% This file is part of Dynare.
%
% Dynare is free software: you can redistribute it and/or modify
% it under the terms of the GNU General Public License as published by
% the Free Software Foundation, either version 3 of the License, or
% (at your option) any later version.
%
% Dynare is distributed in the hope that it will be useful,
% but WITHOUT ANY WARRANTY; without even the implied warranty of
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
% GNU General Public License for more details.
%
% You should have received a copy of the GNU General Public License
% along with Dynare. If not, see <https://www.gnu.org/licenses/>.
C3 = 200;
C4 = 20.2;
C5 = 19.8;
C6 = 180.0;
if nargin==1
fval = zeros(4, 1);
tmp1 = x(2)-x(1)*x(1);
tmp2 = x(4)-x(3)*x(3);
fval(1) = -C3*x(1)*tmp1-(1-x(1));
fval(2) = C3*tmp1+C4*(x(2)-1.0)+C5*(x(4)-1);
fval(3) = -C6*x(3)*tmp2-(1-x(3));
fval(4) = C6*tmp2+C4*(x(4)-1)+C5*(x(2)-1);
if nargout>1
fjac = zeros(4);
tmp1 = x(2)-3*x(1)*x(1);
tmp2 = x(4)-3*x(3)*x(3);
fjac(1,1) = -C3*tmp1+1;
fjac(1,2) = -C3*x(1);
fjac(2,1) = -2*C3*x(1);
fjac(2,2) = C3+C4;
fjac(2,4) = C5;
fjac(3,3) = -C6*tmp2+1;
fjac(3,4) = -C6*x(3);
fjac(4,2) = C5;
fjac(4,3) = -2*C6*x(3);
fjac(4,4) = C6+C4;
end
elseif nargin
if nargout==1
fval = [-3; -1; -3; -1];
else
error('One output is required for initialization mode.')
end
else
error('Wrong number of input arguments.')
end