1function [M2a, M3a, GAMMAa] = amap2_adjust_gamma(M1, M2, M3, GAMMA, weights, method, constraints)
2% Computes a feasible set of characteristics
for the MAP(2) that
is as
3% close as possible to the desired set of characterstics.
5% - M1,M2,M3: moments of the marginal distribution
6% - GAMMA: auto-correlation decay rate
7% - weights: an optional three-element vector with the weights associated
8% to M2, M3, GAMMA. Default
is [10 1 10];
9% - method: the method used to compute the feasible set of characteristics
10% 1) pattern search on M2, M3, GAMMA
11% 2) fit M3 as closely as possible, then perform pattern search
12% on M3, GAMMA (weight of M2
is ignored)
13% 3) (default) prioritize M2, then M3, then GAMMA
14% 4) fit M2 as closely as possible, then perform PSwarm
15% on M3 with bounds [M3_lb(M2), M3_ub(M2)] minimizing the
17% f(x) = weight(2)*(M3a/M3-1)^2 + weight(3)*(GAMMAa/GAMMA-1)^2
18% where GAMMAa = max(GAMMA_LB(M3a), min(GAMMA, GAMMA_UB(M3a))
20% - M2a,M3a: feasible moments of the marginal distribution (M1
is always
22% - GAMMAa: feasible auto-correlation decay rate
24if nargin == 4 || isempty(weights)
25 % ignored with method 3 and 4
38 nonlcon = @nonlcon_safe;
40 nonlcon = @nonlcon_safe2;
42elseif constraints == 2
44 nonlcon = @nonlcon_theoretical;
46 nonlcon = @nonlcon_theoretical2;
49 error('Invalid option');
52if method == 1 || method == 2
54 options = psoptimset(options,'MaxIter', 1e5);
55 options = psoptimset(options,'MaxFunEvals', 1e5);
56 options = psoptimset(options,'PollingOrder', 'Success');
57 options = psoptimset(options,'SearchMethod', { @searchneldermead [] [] });
58 options = psoptimset(options,
'CompleteSearch',
'on');
59 options = psoptimset(options,
'Display',
'iter');
60 options = psoptimset(options,
'TolCon', 1e-100);
63% tolerance
for strict inequalities
68 FEASIBLE = map_scale(amap2_assemble(1,1/3,1/2,2/3,1),M1);
70 FEASIBLE = map_scale(amap2_assemble(1,2/3,1/2,2/3,2),M1);
72 FM2 = map_moment(FEASIBLE,2);
73 FM3 = map_moment(FEASIBLE,3);
74 FGAMMA = map_acf(FEASIBLE,4)/map_acf(FEASIBLE,3);
75 x = patternsearch(@fun,[FM2 FM3 FGAMMA],[],[],[],[],[0 0 -1],[inf inf, 1-tol], nonlcon, options);
80 % force feasibility of the second moment
81 M2a = max(3/2*M1^2, M2);
82 % find a feasible value of the third moment
84 p2 = 3*(n2-2)/(3*n2) * (-2*sqrt(3)/sqrt(12-6*n2) - 1);
85 a2 = (n2-2)/(p2*(1-n2) + sqrt(p2^2 + (2*p2*(n2-2))));
86 l2 = 3*(a2+1)/(a2*p2+1) - (6*a2)/(2+a2*p2*(2*a2+2));
88 if 3/2 <= n2 && n2 < 2
89 FM3 = (l2 + (u2-l2)/2) * M1* M2a;
90 M3_LB = l2 * M1 * M2a;
91 M3_UB = u2 * M1 * M2a;
93 FM3 = (3/2) * M2a^2 / M1 + tol;
97 % null decay rate
is always feasible
99 % find best feasible values of M3 and GAMMA
100 x = patternsearch(@fun2,[FM3 FGAMMA],[],[],[],[],[M3_LB -1],[M3_UB, 1-tol], nonlcon, options);
101 %x = fmincon(@fun2,[FM3 FGAMMA],[],[],[],[],[M3_LB -1],[M3_UB, 1], nonlcon);
105 % priorities are M2 > M3 > GAMMA
106 [M2a,M3a] = aph2_adjust(M1, M2, M3);
107 [lb, ub] = compute_gamma_bounds(M2a, M3a);
108 GAMMAa = max(lb, min(GAMMA, ub));
110 % priorities are M2 > GAMMA > M3
111 % adjust second moment
113 scv = ((M2-M1sq)/M1sq);
115 % (M2 - M1^2)/M1^2 = 1/2
116 % M2 - M1^2 = 1/2 M1^2
119 scva = ((M2a-M1sq)/M1sq);
124 % get bounds on third moment
126 M3_lb = 3*M1^3*(3*scva-1+sqrt(2)*(1-scva)^(3/2));
129 M3_lb = 3/2*M1^3*(1+scva)^2;
132 % compute M3a and GAMMAa
133 if abs(M3_lb - M3_ub) < tol
135 M3a = (M3_lb + M3_ub)/2;
138 % optimization problem formulation
139 problem.Variables = 1;
140 problem.LB = M3_lb + tol;
142 problem.ObjFunction = @gamma_objective;
143 M3a = PSwarm(problem);
144 % compute adjusted gamma
145 [lb, ub] = compute_gamma_bounds(M2a, M3a);
146 GAMMAa = max(lb, min(GAMMA, ub));
149 error(
'Invalid method for adjusting AMAP(2) characteristics');
152 function obj = gamma_objective(M3a)
153 [lb,ub] = compute_gamma_bounds(M2a, M3a);
154 GAMMAa = max(lb, min(GAMMA, ub));
155 obj = weights(2) * (M3a/M3-1)^2 + weights(3) * (GAMMAa/GAMMA-1)^2;
158 function [lb, ub] = compute_gamma_bounds(M2a, M3a)
162 lb = -(N2a*(N3a-6)+6)/(3*N2a-6);
163 ub = -(2* (0.5*(N2a-2)+0.5*sqrt(N2a^2 - 2*N2a*N3a/3))^2 )/(N2a-2);
165 elseif N3a < 9 - 12/N2a
166 lb = -(N2a*(N3a-6)+6)/(3*N2a-6);
169 x1 = sqrt(N2a*(N2a*(18*N2a+N3a*(N3a-18)-27)+24*N3a));
171 lb = (x2-x1+12)/(x2+x1+12);
176 function obj = fun(x)
177 v = (x - [M2 M3 GAMMA]) ./ [M2 M3 GAMMA] .* weights;
181 function obj = fun2(x)
182 v = (x - [M3 GAMMA]) ./ [M3 GAMMA] .* weights(2:3);
186 function [c, ceq] = nonlcon_safe(x)
190 AMAPS = amap2_fit_decay(M1,xM2,xM3,xGAMMA);
199 function [c, ceq] = nonlcon_safe2(x)
202 AMAPS = amap2_fit_decay(M1,M2,xM3,xGAMMA);
211 function [c, ceq] = nonlcon_theoretical(x)
216 % compute normalized moments
220 %fprintf(
'M1 = %f, M2 = %f, M3 = %f, GAMMA = %f\n', M1, xM2, xM3, xGAMMA);
221 %fprintf(
'n2 = %f, n3 = %f\n', n2, n3);
223 % simplifying notations
for moments bounds
224 p2 = 3*(n2-2)/(3*n2) * (-2*sqrt(3)/sqrt(12-6*n2) - 1);
225 a2 = (n2-2)/(p2*(1-n2) + sqrt(p2^2 + (2*p2*(n2-2))));
226 l2 = 3*(a2+1)/(a2*p2+1) - (6*a2)/(2+a2*p2*(2*a2+2));
229 %
if (M2-M1^2)/M1^2 == 1
230 % fprintf(
'Exponential\n');
236 % CONSTRAINT 1: 3/2 <= n
238 if 3/2 <= n2 && n2 < 2
239 % CONSTRAINT 2: l2 <= n3
240 % CONSTRAINT 3: n3 <= u2
244 % CONSTRAINT 2: 3/2 n2 < n3
245 % CONSTRAINT 3: disabled
246 c(2) = 3/2 * n2 - n3 + eps;
250 % BOUNDS FOR AUTOCORRELATION DECAY
251 lb1 = -(n2*(n3-6)+6)/(3*n2-6);
252 ub1 = -(2*(1/2*(n2-2)+1/2*sqrt(n2^2-(2*n2*n3)/3))^2)/(n2-2);
254 tmp2 = sqrt(n2*(n2*(18*n2+n3*(n3-18)-27)+24*n3));
255 lb2 = (tmp1-tmp2+12)/(tmp1+tmp2+12);
256 %fprintf(
'LB1 = %f, LB2 = %f\n', lb1, lb2);
258 %fprintf(
'LB1 <= gamma: %d\n', lb1 <= xGAMMA);
259 %fprintf(
'gamma <= UB: %d\n', xGAMMA <= ub1);
262 elseif n2 > 2 && n3 < 9 - 12/n2
263 %fprintf(
'LB1 <= gamma: %d\n', lb1 <= xGAMMA);
264 %fprintf(
'gamma <= 1: %d\n', xGAMMA <= 1);
267 elseif n2 > 2 && n3 >= 9 - 12/n2
268 %fprintf(
'LB2 <= gamma: %d\n', lb2 <= xGAMMA);
269 %fprintf(
'gamma <= 1: %d\n', xGAMMA <= 1);
278 function [c, ceq] = nonlcon_theoretical2(x)
283 % compute normalized moments
287 % M3
is always feasible because of the bound constraints
289 % BOUNDS FOR AUTOCORRELATION DECAY
290 lb1 = -(n2*(n3-6)+6)/(3*n2-6);
291 ub1 = -(2*(1/2*(n2-2)+1/2*sqrt(n2^2-(2*n2*n3)/3))^2)/(n2-2);
293 tmp2 = sqrt(n2*(n2*(18*n2+n3*(n3-18)-27)+24*n3));
294 lb2 = (tmp1-tmp2+12)/(tmp1+tmp2+12);
295 %fprintf(
'LB1 = %f, LB2 = %f\n', lb1, lb2);
297 %fprintf(
'LB1 <= gamma: %d\n', lb1 <= xGAMMA);
298 %fprintf(
'gamma <= UB: %d\n', xGAMMA <= ub1);
301 elseif n2 > 2 && n3 < 9 - 12/n2
302 %fprintf(
'LB1 <= gamma: %d\n', lb1 <= xGAMMA);
303 %fprintf(
'gamma <= 1: %d\n', xGAMMA <= 1);
306 elseif n2 > 2 && n3 >= 9 - 12/n2
307 %fprintf(
'LB2 <= gamma: %d\n', lb2 <= xGAMMA);
308 %fprintf(
'gamma <= 1: %d\n', xGAMMA <= 1);