Image Interpolation with Contour Stencils
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros
cwinterp.c
Go to the documentation of this file.
1 
60 #include <stdio.h>
61 #include <stdlib.h>
62 #include <string.h>
63 #include <math.h>
64 
65 #include "basic.h"
66 #include "cwinterp.h"
67 #include "fitsten.h"
68 #include "invmat.h"
69 #include "nninterp.h"
70 #include "drawline.h"
71 
72 
74 #define NUMSTENCILS 8
75 
77 #define NEIGHRADIUS 1
78 #define NEIGHDIAMETER (2*NEIGHRADIUS+1)
79 #define NUMNEIGH (NEIGHDIAMETER*NEIGHDIAMETER)
80 
89 #define CORRECTION_IGNOREBITS 3
90 
92 #define INPUT_FRACBITS 8
93 
94 #define PSI_FRACBITS 12
95 
96 #define OUTPUT_FRACBITS (INPUT_FRACBITS + PSI_FRACBITS)
97 
99 #define FIXED_ONE(N) (1 << (N))
100 
101 #define FIXED_HALF(N) (1 << ((N) - 1))
102 
103 
114 #define PIXEL_STRIDE 3
115 
116 
117 /* Generic macros */
118 
120 #define CLAMP(X,A,B) (((X) < (A)) ? (A) : (((X) > (B)) ? (B) : (X)))
121 
123 #define ROUNDCLAMP(X,A,B) (((X) < (A)) ? (A) : (((X) > (B)) ? (B) : ROUND(X)))
124 
125 #ifndef M_2PI
126 
127 #define M_2PI 6.283185307179586476925286766559
128 #endif
129 
130 
131 /* Orientation in radians of each stencil */
132 static const double StencilOrientation[NUMSTENCILS] = {-1.178097245096172,
133  -0.785398163397448, -0.392699081698724, 0.0, 0.392699081698724,
134  0.785398163397448, 1.178097245096172, 1.570796326794897};
135 
136 
138 static double Psf(double x, double y, cwparams Param)
139 {
140  double SigmaSqr = Param.PsfSigma*Param.PsfSigma;
141  return exp(-(x*x + y*y)/(2.0*SigmaSqr))/(M_2PI*SigmaSqr);
142 }
143 
144 
146 static double Phi(double x, double y,
147  double theta, double PhiSigmaTangent, double PhiSigmaNormal)
148 {
149  double t, n;
150 
151  /* Oriented Gaussian */
152  t = (cos(theta)*x + sin(theta)*y) / PhiSigmaTangent;
153  n = (-sin(theta)*x + cos(theta)*y) / PhiSigmaNormal;
154 
155  return exp(-0.5*(t*t + n*n));
156 }
157 
158 
160 static float CubicBSpline(float x)
161 {
162  x = (float)fabs(x);
163 
164  if(x < 1)
165  return (x/2 - 1)*x*x + 0.66666666666666667f;
166  else if(x < 2)
167  {
168  x = 2 - x;
169  return x*x*x/6;
170  }
171  else
172  return 0;
173 }
174 
175 
177 static double Window(double x, double y)
178 {
179  double Temp;
180 
181  x *= 2.0/(NEIGHRADIUS + 1.0);
182  y *= 2.0/(NEIGHRADIUS + 1.0);
183 
184  /* Cubic B-spline */
185  if(-2.0 < x && x < 2.0 && -2.0 < y && y < 2.0)
186  {
187  x = fabs(x);
188  Temp = fabs(1.0 - x);
189  x = 1.0 - x + (x*x*x - 2.0*Temp*Temp*Temp)/6.0;
190 
191  y = fabs(y);
192  Temp = fabs(1.0 - y);
193  y = 1.0 - y + (y*y*y - 2.0*Temp*Temp*Temp)/6.0;
194  return (x * y) * 4.0/((NEIGHRADIUS + 1.0)*(NEIGHRADIUS + 1.0));
195  }
196  else
197  return 0.0;
198 }
199 
200 
202 static void QuadraturePoint(double *Weight, double *Abscissa, int Index, int NumPanels)
203 {
204  switch(Index % 3)
205  {
206  case 0:
207  *Weight = (Index == 0 || Index == NumPanels)? 0.25 : 0.5;
208  *Abscissa = Index;
209  break;
210  case 1:
211  *Weight = 1.25;
212  /* Abscissa location is Index - (3.0/sqrt(5.0) - 1.0)/2.0 */
213  *Abscissa = Index - 1.70820393249936919e-1;
214  break;
215  case 2:
216  *Weight = 1.25;
217  /* Abscissa location is Index + (3.0/sqrt(5.0) - 1.0)/2.0 */
218  *Abscissa = Index + 1.70820393249936919e-1;
219  break;
220  }
221 }
222 
223 
225 static double PsfPhiConvolution(int x, int y, double Theta, cwparams Param)
226 {
227  /* Integrate over the square [-R,R]x[-R,R] */
228  const double R = 4.0*Param.PsfSigma;
229  /* Number of panels along each dimension, must be divisible by 3 */
230  const int NumPanels = 3*16;
231  const double PanelSize = 2.0*R/NumPanels;
232  double Integral = 0.0, Slice = 0.0;
233  double u, v, wu, wv;
234  int IndexX, IndexY;
235 
236 
237  /* Specially handle the case where PSF is Dirac delta */
238  if(Param.PsfSigma == 0.0)
239  return Phi(x, y, Theta, Param.PhiSigmaTangent, Param.PhiSigmaNormal);
240 
241  /* Approximate 2D integral */
242  for(IndexY = 0; IndexY <= NumPanels; IndexY++)
243  {
244  QuadraturePoint(&wv, &v, IndexY, NumPanels);
245  v = PanelSize*v - R;
246 
247  for(Slice = 0.0, IndexX = 0; IndexX <= NumPanels; IndexX++)
248  {
249  QuadraturePoint(&wu, &u, IndexX, NumPanels);
250  u = PanelSize*u - R;
251  Slice += wu*( Psf(u, v, Param) *
252  Phi(x - u, y - v, Theta, Param.PhiSigmaTangent,
253  Param.PhiSigmaNormal) );
254  }
255 
256  Integral += wv*Slice;
257  }
258 
259  Integral *= PanelSize*PanelSize;
260  return Integral;
261 }
262 
263 
265 static int ComputeMatrices(double *InverseA, cwparams Param)
266 {
267  double *A = NULL;
268  int m, mx, my, n, nx, ny, S;
269  int Status = 0;
270 
271 
272  if(!(A = (double *)Malloc(sizeof(double)*NUMNEIGH*NUMNEIGH)))
273  goto Catch;
274 
275  for(S = 0; S < NUMSTENCILS; S++)
276  {
277  for(ny = -NEIGHRADIUS, n = 0; ny <= NEIGHRADIUS; ny++)
278  for(nx = -NEIGHRADIUS; nx <= NEIGHRADIUS; nx++, n++)
279  for(my = -NEIGHRADIUS, m = 0; my <= NEIGHRADIUS; my++)
280  for(mx = -NEIGHRADIUS; mx <= NEIGHRADIUS; mx++, m++)
281  {
282  A[m + NUMNEIGH*n] = PsfPhiConvolution(mx - nx, my - ny,
283  StencilOrientation[S], Param);
284  }
285 
286  /* Compute the inverse of A */
287  if(!InvertMatrix(InverseA + S*(NUMNEIGH*NUMNEIGH), A, NUMNEIGH))
288  goto Catch;
289  }
290 
291  Status = 1;
292 Catch: /* This label is used for error handling. If something went wrong
293  above (which may be out of memory or a computation error), then
294  execution jumps to this point to clean up and exit. */
295  Free(A);
296  return Status;
297 }
298 
299 
300 #include "imageio.h"
301 
319 int32_t *PreCWInterp(cwparams Param)
320 {
321  int32_t *Psi = NULL;
322  const int ScaleFactor = (int)ceil(Param.ScaleFactor);
323  const int SupportRadius = (NEIGHRADIUS+1)*ScaleFactor - 1;
324  const int SupportWidth = 2*SupportRadius + 1;
325  const int SupportSize = SupportWidth*SupportWidth;
326  double *InverseA = NULL;
327  double x, y, Wxy, Psi0, Psim, XStart, YStart, WSum;
328  int S, sx, sy, i, m0, m, mx, my, n, nx, ny, Success = 0;
329 
330 
331  if(!(Psi = (int32_t *)Malloc(sizeof(int32_t)*SupportSize*NUMNEIGH*NUMSTENCILS))
332  || !(InverseA = (double *)Malloc(sizeof(double)*NUMNEIGH*NUMNEIGH*NUMSTENCILS)))
333  goto Catch;
334 
335  /* Compute the matrices, the results are stored in InverseA. */
336  if(!ComputeMatrices(InverseA, Param))
337  goto Catch;
338 
339  if(Param.CenteredGrid)
340  {
341  XStart = (1/Param.ScaleFactor - 1)/2;
342  YStart = (1/Param.ScaleFactor - 1)/2;
343  }
344  else
345  XStart = YStart = 0;
346 
348 
349  /* Precompute the samples of the Psi functions */
350  for(S = 0; S < NUMSTENCILS; S++)
351  for(i = 0, sy = -SupportRadius; sy <= SupportRadius; sy++)
352  for(sx = -SupportRadius; sx <= SupportRadius; sx++, i++)
353  {
354  /* Compute the sum of window translates. This sum should be
355  exactly constant, but there can be small variations. We
356  divide the Psi samples computed below by WSum to compensate.
357  */
358  for(ny = -(int)floor((sy + SupportRadius)/ScaleFactor), WSum = 0;
359  ny <= (2*NEIGHRADIUS + 1) && sy + ny*ScaleFactor <= SupportRadius; ny++)
360  for(nx = -(int)floor((sx + SupportRadius)/ScaleFactor);
361  nx <= (2*NEIGHRADIUS + 1) && sx + nx*ScaleFactor <= SupportRadius; nx++)
362  {
363  x = XStart + nx + ((double)sx)/((double)ScaleFactor);
364  y = YStart + ny + ((double)sy)/((double)ScaleFactor);
365  WSum += ROUND(Window(x,y)*FIXED_ONE(PSI_FRACBITS))
366  / (double)FIXED_ONE(PSI_FRACBITS);
367  }
368 
369  x = XStart + ((double)sx)/((double)ScaleFactor);
370  y = YStart + ((double)sy)/((double)ScaleFactor);
371  Psi0 = Wxy = Window(x, y);
372 
373  for(my = -NEIGHRADIUS, m = 0; my <= NEIGHRADIUS; my++)
374  for(mx = -NEIGHRADIUS; mx <= NEIGHRADIUS; mx++, m++)
375  {
376  if(m != m0)
377  {
378  Psim = 0.0;
379 
380  for(ny = -NEIGHRADIUS, n = 0; ny <= NEIGHRADIUS; ny++)
381  for(nx = -NEIGHRADIUS; nx <= NEIGHRADIUS; nx++, n++)
382  Psim += InverseA[m + NUMNEIGH*(n + NUMNEIGH*S)]
383  * Phi(x - nx, y - ny,
384  StencilOrientation[S], Param.PhiSigmaTangent,
385  Param.PhiSigmaNormal);
386 
387  Psim *= Wxy;
388  Psi0 -= Psim;
389 
390  Psi[i + SupportSize*(m + NUMNEIGH*S)] =
391  (int32_t)ROUND((Psim/WSum)*FIXED_ONE(PSI_FRACBITS));
392  }
393  }
394 
395  Psi[i + SupportSize*(m0 + NUMNEIGH*S)] =
396  (int32_t)ROUND((Psi0/WSum)*FIXED_ONE(PSI_FRACBITS));
397  }
398 
399  Success = 1;
400 Catch: /* This label is used for error handling. If something went wrong
401  above (which may be out of memory or a computation error), then
402  execution jumps to this point to clean up and exit. */
403  Free(InverseA);
404  if(!Success && Psi)
405  {
406  Free(Psi);
407  Psi = NULL;
408  }
409  return Psi;
410 }
411 
412 
428 static void CWFirstPass(int32_t *Interpolation, int ScaleFactor, const int32_t *Input,
429  int InputWidth, int InputHeight, const int *Stencil, const int32_t *Psi)
430 {
431  const int32_t *PsiPtr, *SrcWindow;
432  int32_t *DestWindow;
433 
434  const int Pad = 2;
435  const int SampleRange = (NEIGHRADIUS+1)*ScaleFactor - 1;
436  const int SampleWidth = 2*SampleRange + 1;
437 
438  const int OutputWidth = ScaleFactor*InputWidth;
439  const int DestWindowJump = PIXEL_STRIDE*(OutputWidth - SampleWidth);
440  const int DestStep = PIXEL_STRIDE*ScaleFactor;
441  const int DestJump = PIXEL_STRIDE*(ScaleFactor-1)*OutputWidth + 2*Pad*DestStep;
442  const int SrcWindowJump = PIXEL_STRIDE*(InputWidth - NEIGHDIAMETER);
443  const int SrcJump = 2*PIXEL_STRIDE*Pad;
444  const int StencilJump = 2*Pad;
445 
446  int x, y, NeighX, NeighY, SampleX, SampleY;
447  int32_t cr, cg, cb;
448 
449 
450  Interpolation += PIXEL_STRIDE*(Pad*ScaleFactor - SampleRange)*(1 + OutputWidth);
451  Input += PIXEL_STRIDE*(Pad - NEIGHRADIUS)*(1 + InputWidth);
452  Stencil += Pad*(1 + InputWidth);
453 
454  for(y = InputHeight - 2*Pad; y; y--,
455  Stencil += StencilJump, Input += SrcJump, Interpolation += DestJump)
456  for(x = InputWidth - 2*Pad; x; x--,
457  Stencil++, Input += PIXEL_STRIDE, Interpolation += DestStep)
458  {
459  PsiPtr = Psi + *Stencil;
460  SrcWindow = Input;
461 
462  for(NeighY = NEIGHDIAMETER; NeighY; NeighY--, SrcWindow += SrcWindowJump)
463  for(NeighX = NEIGHDIAMETER; NeighX; NeighX--, SrcWindow += PIXEL_STRIDE)
464  {
465  cr = SrcWindow[0];
466  cg = SrcWindow[1];
467  cb = SrcWindow[2];
468  DestWindow = Interpolation;
469  SampleY = SampleWidth;
470 
471  for(SampleY = SampleWidth; SampleY;
472  SampleY--, DestWindow += DestWindowJump, PsiPtr += SampleWidth)
473  for(SampleX = 0; SampleX < SampleWidth;
474  SampleX++, DestWindow += PIXEL_STRIDE)
475  {
476  int32_t Temp = PsiPtr[SampleX];
477  DestWindow[0] += cr * Temp;
478  DestWindow[1] += cg * Temp;
479  DestWindow[2] += cb * Temp;
480  }
481  }
482  }
483 }
484 
485 
501 static void CWRefinementPass(int32_t *Interpolation, int ScaleFactor,
502  const int32_t *Residual, int InputWidth, int InputHeight,
503  const int *Stencil, const int32_t *Sample)
504 {
505  const int Pad = 4;
506 
507  const int32_t *SamplePtr, *SrcWindow;
508  int32_t *DestWindow;
509  const int SampleRange = (NEIGHRADIUS+1)*ScaleFactor - 1;
510  const int SampleWidth = 2*SampleRange + 1;
511  const int SampleSize = SampleWidth*SampleWidth;
512 
513  const int OutputWidth = ScaleFactor*InputWidth;
514  const int DestWindowJump = PIXEL_STRIDE*(OutputWidth - SampleWidth);
515  const int DestStep = PIXEL_STRIDE*ScaleFactor;
516  const int DestJump = PIXEL_STRIDE*(ScaleFactor-1)*OutputWidth + 2*Pad*DestStep;
517  const int SrcWindowJump = PIXEL_STRIDE*(InputWidth - NEIGHDIAMETER);
518  const int SrcJump = 2*PIXEL_STRIDE*Pad;
519  const int StencilJump = 2*Pad;
520 
521  int x, y, NeighX, NeighY, SampleX, SampleY;
522  int32_t cr, cg, cb;
523 
524 
525  Interpolation += PIXEL_STRIDE*(Pad*ScaleFactor - SampleRange)*(1 + OutputWidth);
526  Residual += PIXEL_STRIDE*(Pad - NEIGHRADIUS)*(1 + InputWidth);
527  Stencil += Pad*(1 + InputWidth);
528 
529  for(y = InputHeight - 2*Pad; y; y--,
530  Stencil += StencilJump, Residual += SrcJump, Interpolation += DestJump)
531  for(x = InputWidth - 2*Pad; x; x--,
532  Stencil++, Residual += PIXEL_STRIDE, Interpolation += DestStep)
533  {
534  SamplePtr = Sample + *Stencil;
535  SrcWindow = Residual;
536 
537  for(NeighY = NEIGHDIAMETER; NeighY; NeighY--, SrcWindow += SrcWindowJump)
538  for(NeighX = NEIGHDIAMETER; NeighX; NeighX--, SrcWindow += PIXEL_STRIDE)
539  {
540  cr = SrcWindow[0];
541  cg = SrcWindow[1];
542  cb = SrcWindow[2];
543 
544  if( ((cr >> CORRECTION_IGNOREBITS) && (-cr >> CORRECTION_IGNOREBITS))
545  || ((cg >> CORRECTION_IGNOREBITS) && (-cg >> CORRECTION_IGNOREBITS))
546  || ((cb >> CORRECTION_IGNOREBITS) && (-cb >> CORRECTION_IGNOREBITS)) )
547  {
548  DestWindow = Interpolation;
549  SampleY = SampleWidth;
550 
551  for(SampleY = SampleWidth; SampleY;
552  SampleY--, DestWindow += DestWindowJump, SamplePtr += SampleWidth)
553  for(SampleX = 0; SampleX < SampleWidth;
554  SampleX++, DestWindow += PIXEL_STRIDE)
555  {
556  int32_t Temp = SamplePtr[SampleX];
557  DestWindow[0] += cr * Temp;
558  DestWindow[1] += cg * Temp;
559  DestWindow[2] += cb * Temp;
560  }
561  }
562  else
563  SamplePtr += SampleSize;
564  }
565  }
566 }
567 
568 
575 static int ConstExtension(int N, int i)
576 {
577  if(i < 0)
578  return 0;
579  else if(i >= N)
580  return N - 1;
581  else
582  return i;
583 }
584 
585 
586 static float Sqr(float x)
587 {
588  return x*x;
589 }
590 
591 
593 static int32_t CWResidual(int32_t *Residual, const int32_t *Interpolation,
594  const int32_t *Input, int CoarseWidth, int CoarseHeight, cwparams Param)
595 {
596  int Pad = 4;
597  const int ScaleFactor = (int)ceil(Param.ScaleFactor);
598  const int InterpWidth = ScaleFactor*CoarseWidth;
599  const int InterpHeight = ScaleFactor*CoarseHeight;
600  const int CoarseStride = 3*CoarseWidth;
601  const float PsfRadius = (float)(4*Param.PsfSigma*ScaleFactor);
602  const int PsfWidth = (int)ceil(2*PsfRadius);
603  float *Temp = NULL, *PsfBuf = NULL;
604  float ExpDenom, Weight, Sum[3], DenomSum;
605  float XStart, YStart, X, Y;
606  int IndexX0, IndexY0, SrcOffset, DestOffset;
607  int x, y, i, n, c, Success = 0;
608  int x1, x2;
609  int32_t ResNorm = 0;
610 
611 
612  if(!(Temp = (float *)Malloc(sizeof(float)*3*CoarseWidth*InterpHeight))
613  || !(PsfBuf = (float *)Malloc(sizeof(float)*PsfWidth)))
614  goto Catch;
615 
616  if(Param.CenteredGrid)
617  {
618  XStart = (1.0f/ScaleFactor - 1)/2;
619  YStart = (1.0f/ScaleFactor - 1)/2;
620  }
621  else
622  XStart = YStart = 0;
623 
624  if(Param.PsfSigma)
625  ExpDenom = 2 * Sqr((float)(Param.PsfSigma*ScaleFactor));
626  else
627  ExpDenom = 2 * Sqr(1e-2f*ScaleFactor);
628 
629  if(Pad < ScaleFactor)
630  Pad = ScaleFactor;
631 
632  for(x = 0; x < CoarseWidth; x++)
633  {
634  X = (-XStart + x)*ScaleFactor;
635  IndexX0 = (int)ceil(X - PsfRadius);
636 
637  /* Evaluate the PSF */
638  for(n = 0; n < PsfWidth; n++)
639  PsfBuf[n] = (float)exp(-Sqr(X - (IndexX0 + n)) / ExpDenom);
640 
641  for(y = 0, SrcOffset = 0, DestOffset = 3*x; y < InterpHeight;
642  y++, SrcOffset += InterpWidth, DestOffset += CoarseStride)
643  {
644  Sum[0] = Sum[1] = Sum[2] = DenomSum = 0;
645 
646  for(n = 0; n < PsfWidth; n++)
647  {
648  Weight = PsfBuf[n];
649  DenomSum += Weight;
650  i = 3*(ConstExtension(InterpWidth, IndexX0 + n) + SrcOffset);
651 
652  for(c = 0; c < 3; c++)
653  Sum[c] += Weight * Interpolation[i + c];
654  }
655 
656  for(c = 0; c < 3; c++)
657  Temp[DestOffset + c] = Sum[c] / DenomSum;
658  }
659  }
660 
661  x1 = 3*Pad;
662  x2 = CoarseStride - 3*Pad;
663 
664  for(y = 0; y < CoarseHeight; y++,
665  Residual += CoarseStride, Input += CoarseStride)
666  {
667  if(!(y >= Pad && y < CoarseHeight-Pad))
668  continue;
669 
670  Y = (-YStart + y)*ScaleFactor;
671  IndexY0 = (int)ceil(Y - PsfRadius);
672 
673  /* Evaluate the PSF */
674  for(n = 0; n < PsfWidth; n++)
675  PsfBuf[n] = (float)exp(-Sqr(Y - (IndexY0 + n)) / ExpDenom);
676 
677  for(x = x1; x < x2; x += 3)
678  {
679  Sum[0] = Sum[1] = Sum[2] = DenomSum = 0;
680 
681  for(n = 0; n < PsfWidth; n++)
682  {
683  SrcOffset = x + CoarseStride*ConstExtension(InterpHeight, IndexY0 + n);
684  Weight = PsfBuf[n];
685  DenomSum += Weight;
686 
687  for(c = 0; c < 3; c++)
688  Sum[c] += Weight * Temp[SrcOffset + c];
689  }
690 
691  DenomSum *= FIXED_ONE(PSI_FRACBITS);
692 
693  for(c = 0; c < 3; c++)
694  {
695  Sum[c] = Input[x + c] - Sum[c] / DenomSum;
696  Residual[x + c] = (int32_t)ROUND(Sum[c]);
697 
698  if(abs(Residual[x + c]) > ResNorm)
699  ResNorm = abs(Residual[x + c]);
700  }
701  }
702  }
703 
704  Success = 1;
705 Catch:
706  Free(PsfBuf);
707  Free(Temp);
708  return (Success) ? ResNorm : -1;
709 }
710 
711 
730 static void ConvertInput(int32_t *FixedRgb, const uint32_t *Input, int InputWidth,
731  int InputHeight, int Padding)
732 {
733  const uint8_t *InputPtr = (uint8_t *)Input;
734  const int InputStride = 4*InputWidth;
735  const int Stride = PIXEL_STRIDE*(InputWidth + 2*Padding);
736  int32_t r, g, b;
737  int i, Row;
738 
739 
740  FixedRgb += Padding*Stride;
741 
742  for(Row = InputHeight; Row; Row--, InputPtr += InputStride)
743  {
744  r = ((int32_t)InputPtr[0]) << INPUT_FRACBITS;
745  g = ((int32_t)InputPtr[1]) << INPUT_FRACBITS;
746  b = ((int32_t)InputPtr[2]) << INPUT_FRACBITS;
747 
748  /* Pad left side by copying pixel */
749  for(i = Padding; i; i--)
750  {
751  *(FixedRgb++) = r;
752  *(FixedRgb++) = g;
753  *(FixedRgb++) = b;
754  }
755 
756  /* Convert the interior of the image */
757  for(i = 0; i < InputStride; i += 4)
758  {
759  *(FixedRgb++) = ((int32_t)InputPtr[i+0]) << INPUT_FRACBITS;
760  *(FixedRgb++) = ((int32_t)InputPtr[i+1]) << INPUT_FRACBITS;
761  *(FixedRgb++) = ((int32_t)InputPtr[i+2]) << INPUT_FRACBITS;
762  }
763 
764  r = ((int32_t)InputPtr[i-4]) << INPUT_FRACBITS;
765  g = ((int32_t)InputPtr[i-3]) << INPUT_FRACBITS;
766  b = ((int32_t)InputPtr[i-2]) << INPUT_FRACBITS;
767 
768  /* Pad right side by copying pixel */
769  for(i = Padding; i; i--)
770  {
771  *(FixedRgb++) = r;
772  *(FixedRgb++) = g;
773  *(FixedRgb++) = b;
774  }
775  }
776 
777  /* Pad bottom by copying rows */
778  for(Row = Padding; Row; Row--, FixedRgb += Stride)
779  memcpy(FixedRgb, FixedRgb - Stride, sizeof(int32_t)*Stride);
780 
781  FixedRgb -= Stride*(InputHeight + Padding);
782 
783  /* Pad top by coping rows */
784  for(Row = Padding; Row; Row--, FixedRgb -= Stride)
785  memcpy(FixedRgb - Stride, FixedRgb, sizeof(int32_t)*Stride);
786 }
787 
788 
812 static void ConvertOutput(uint32_t *Output, int OutputWidth, int OutputHeight,
813  const int32_t *FixedRgb, int Width)
814 {
815  uint8_t *OutputPtr = (uint8_t *)Output;
816  const int CroppedStride = PIXEL_STRIDE*OutputWidth, Stride = PIXEL_STRIDE*Width;
817  int32_t r, g, b;
818  int i, Row;
819 
820 
821  for(Row = OutputHeight; Row; Row--, FixedRgb += Stride)
822  for(i = 0; i < CroppedStride; i += PIXEL_STRIDE)
823  {
824  /* Convert fixed-point values to integer */
825  r = (FixedRgb[i+0] + FIXED_HALF(OUTPUT_FRACBITS)) >> OUTPUT_FRACBITS;
826  g = (FixedRgb[i+1] + FIXED_HALF(OUTPUT_FRACBITS)) >> OUTPUT_FRACBITS;
827  b = (FixedRgb[i+2] + FIXED_HALF(OUTPUT_FRACBITS)) >> OUTPUT_FRACBITS;
828 
829  /* Clamp range to [0, 255] and store in Output */
830  *(OutputPtr++) = CLAMP(r, 0, 255);
831  *(OutputPtr++) = CLAMP(g, 0, 255);
832  *(OutputPtr++) = CLAMP(b, 0, 255);
833  *(OutputPtr++) = 0xFF;
834  }
835 }
836 
837 
847 int CWInterp(uint32_t *Output, const uint32_t *Input,
848  int InputWidth, int InputHeight, const int32_t *Psi, cwparams Param)
849 {
850  const int ScaleFactor = (int)ceil(Param.ScaleFactor);
851  const int SupportRadius = (NEIGHRADIUS+1)*ScaleFactor - 1;
852  const int SupportWidth = 2*SupportRadius + 1;
853  const int SupportSize = SupportWidth*SupportWidth;
854  const int StencilMul = NUMNEIGH*SupportSize;
855  int *Stencil = NULL;
856  int32_t *InputFixed = NULL, *OutputFixed = NULL, *Residual = NULL;
857  unsigned long StartTime, StopTime;
858  int i, PadInput, pw, ph, Success = 0;
859  int32_t ResNorm;
860 
861 
862  /* Iterative refinement is unnecessary if PSF is the Dirac delta */
863  if(Param.PsfSigma == 0.0)
864  Param.RefinementSteps = 0;
865 
866  PadInput = 4 + (ScaleFactor + 1)/2;
867  pw = InputWidth + 2*PadInput;
868  ph = InputHeight + 2*PadInput;
869 
870  if( !(OutputFixed = (int32_t *)Malloc(sizeof(int32_t)*
871  PIXEL_STRIDE*pw*ScaleFactor*ph*ScaleFactor))
872  || !(InputFixed = (int32_t *)Malloc(sizeof(int32_t)*PIXEL_STRIDE*pw*ph))
873  || !(Residual = (int32_t *)Malloc(sizeof(int32_t)*PIXEL_STRIDE*pw*ph))
874  || !(Stencil = (int *)Malloc(sizeof(int)*pw*ph)) )
875  goto Catch;
876 
877  /* Start timing */
878  StartTime = Clock();
879 
880  /* Convert 32-bit RGBA pixels to integer array */
881  ConvertInput(InputFixed, Input, InputWidth, InputHeight, PadInput);
882 
883  /* Select the best-fitting contour stencils */
884  if(!FitStencils(Stencil, InputFixed, pw, ph, StencilMul))
885  goto Catch;
886 
887  memset(OutputFixed, 0, sizeof(int32_t)*
888  3*pw*ScaleFactor*ph*ScaleFactor);
889  memset(Residual, 0, sizeof(int32_t)*3*pw*ph);
890 
891  printf("\n Iteration Residual norm\n -------------------------\n");
892 
893  /* First interpolation pass */
894  CWFirstPass(OutputFixed, ScaleFactor, InputFixed, pw, ph, Stencil, Psi);
895 
896  /* Iterative refinement */
897  for(i = 1; i <= Param.RefinementSteps; i++)
898  {
899  /* Compute the residual */
900  if((ResNorm = CWResidual(Residual, OutputFixed, InputFixed,
901  pw, ph, Param)) < 0.0)
902  goto Catch;
903 
904  printf(" %8d %15.8f\n", i, ResNorm/(255.0*256.0));
905 
906  /* Interpolation refinement pass */
907  CWRefinementPass(OutputFixed, ScaleFactor, Residual, pw, ph, Stencil, Psi);
908  }
909 
910  /* Convert output integer array to 32-bit RGBA */
911  ConvertOutput(Output, InputWidth*ScaleFactor, InputHeight*ScaleFactor,
912  OutputFixed + PIXEL_STRIDE*PadInput*ScaleFactor*(1 + pw*ScaleFactor),
913  pw*ScaleFactor);
914 
915  /* The final interpolation is now complete, stop timing. */
916  StopTime = Clock();
917 
918  /* Compute the residual norm of the final interpolation. This
919  computation is not included in the CPU timing since it is for
920  information purposes only. */
921  ResNorm = CWResidual(Residual, OutputFixed, InputFixed, pw, ph, Param);
922  printf(" %8d %15.8f\n\n", Param.RefinementSteps + 1, ResNorm/(255.0*256.0));
923 
924  /* Display the CPU time spent performing the interpolation. */
925  printf(" CPU time: %.3f s\n\n", 0.001*(StopTime - StartTime));
926 
927  Success = 1;
928 
929 Catch: /* This label is used for error handling. If something went wrong
930  above (which may be out of memory or a computation error), then
931  execution jumps to this point to clean up and exit. */
932  Free(Stencil);
933  Free(Residual);
934  Free(InputFixed);
935  Free(OutputFixed);
936  return Success;
937 }
938 
939 
941 #define ROUNDCLAMP(X,A,B) (((X) < (A)) ? (A) : (((X) > (B)) ? (B) : ROUND(X)))
942 
943 
944 /* The following parameters define the number of fractional bits used for
945  signed 32-bit fixedpoint arithmetic in CWSynth2Fixed. These parameters
946  should be large enough for reasonable precision but small enough to
947  avoid overflow. Additionally, the implementation constraints on
948  choosing these parameters are
949 
950  WINDOW_FRACBITS >= UK_FRACBITS,
951  TRIG_FRACBITS + XY_FRACBITS - PHITN_FRACBITS >= 1,
952  COEFF_FRACBITS + PHI_FRACBITS - UK_FRACBITS >= 1.
953  */
954 #define XY_FRACBITS 15
955 #define TRIG_FRACBITS 10
956 #define PHITN_FRACBITS 9
957 #define PHI_FRACBITS 10
958 #define COEFF_FRACBITS 10
959 #define WINDOW_FRACBITS 10
960 #define UK_FRACBITS 10
961 
962 #define ROUND_FIXED(X,N) (((X) + FIXED_HALF(N)) >> (N))
963 #define FLOAT_TO_FIXED(X,N) ((int32_t)ROUND((X) * FIXED_ONE(N)))
964 
966 int CWArbitraryInterp(uint32_t *Output, int OutputWidth, int OutputHeight,
967  const int32_t *Input, int InputWidth, int InputHeight,
968  const int *Stencil, const double *InverseA, cwparams Param)
969 {
970  /*int (*Extension)(int, int) = ExtensionMethod[Param.Boundary];*/
971  int (*Extension)(int, int) = ConstExtension;
972  const int InputNumEl = 3*InputWidth*InputHeight;
973  const int ExpTableSize = 1024;
974  const double ExpArgScale = 37.0236;
975  const double PhiTScale = sqrt(ExpArgScale/2)/Param.PhiSigmaTangent;
976  const double PhiNScale = sqrt(ExpArgScale/2)/Param.PhiSigmaNormal;
977 
978  int32_t *Coeff = NULL, *CoeffPtr, *ExpTable = NULL;
979 
980  float X, Y, XStart, YStart;
981  float Temp, cr[NUMNEIGH], cg[NUMNEIGH], cb[NUMNEIGH];
982  int32_t v0[3], v[3], WindowWeight, WindowWeightX[4], WindowWeightY[4];
983  int32_t Xpf, Ypf, Weight, uk[3], u[3], DenomSum;
984  int32_t CosTableTf[NUMSTENCILS], SinTableTf[NUMSTENCILS];
985  int32_t CosTableNf[NUMSTENCILS], SinTableNf[NUMSTENCILS];
986  int32_t Pixel;
987  int i, k, x, y, m, n, mx, my, nx, ny, S, Success = 0;
988  int ix, iy, Cur, Offset;
989 
990 
991  if(!(Coeff = (int32_t *)Malloc(sizeof(int32_t)*(NUMNEIGH + 1)*InputNumEl))
992  || !(ExpTable = (int32_t *)Malloc(sizeof(int32_t)*ExpTableSize)))
993  goto Catch;
994 
995  if(Param.CenteredGrid)
996  {
997  XStart = (float)(1/Param.ScaleFactor - 1)/2;
998  YStart = (float)(1/Param.ScaleFactor - 1)/2;
999  }
1000  else
1001  XStart = YStart = 0;
1002 
1003  for(S = 0; S < NUMSTENCILS; S++)
1004  {
1005  CosTableTf[S] = FLOAT_TO_FIXED(
1006  PhiTScale * cos(StencilOrientation[S]), TRIG_FRACBITS);
1007  SinTableTf[S] = FLOAT_TO_FIXED(
1008  PhiTScale * sin(StencilOrientation[S]), TRIG_FRACBITS);
1009  CosTableNf[S] = FLOAT_TO_FIXED(
1010  PhiNScale * cos(StencilOrientation[S]), TRIG_FRACBITS);
1011  SinTableNf[S] = FLOAT_TO_FIXED(
1012  PhiNScale * sin(StencilOrientation[S]), TRIG_FRACBITS);
1013  }
1014 
1015  for(i = 0; i < ExpTableSize; i++)
1016  ExpTable[i] = FLOAT_TO_FIXED(exp( -(double)(i + 0.5f)/ExpArgScale), PHI_FRACBITS);
1017 
1018  for(y = 0, k = 0, CoeffPtr = Coeff; y < InputHeight; y++)
1019  for(x = 0; x < InputWidth; x++, k++, CoeffPtr += 3*(NUMNEIGH + 1))
1020  {
1021  S = NUMNEIGH*Stencil[k];
1022 
1023  v0[0] = Input[3*k + 0];
1024  v0[1] = Input[3*k + 1];
1025  v0[2] = Input[3*k + 2];
1026 
1027  for(m = 0; m < NUMNEIGH; m++)
1028  cr[m] = 0;
1029  for(m = 0; m < NUMNEIGH; m++)
1030  cg[m] = 0;
1031  for(m = 0; m < NUMNEIGH; m++)
1032  cb[m] = 0;
1033 
1034  for(ny = -NEIGHRADIUS, n = 0; ny <= NEIGHRADIUS; ny++)
1035  {
1036  Offset = InputWidth*Extension(InputHeight, y + ny);
1037 
1038  for(nx = -NEIGHRADIUS; nx <= NEIGHRADIUS; nx++, n++)
1039  {
1040  if(n != 1 + (2*NEIGHRADIUS + 1))
1041  {
1042  i = 3*(Extension(InputWidth, x + nx) + Offset);
1043  v[0] = Input[i + 0];
1044  v[1] = Input[i + 1];
1045  v[2] = Input[i + 2];
1046 
1047  for(m = 0; m < NUMNEIGH; m++)
1048  {
1049  Temp = (float)InverseA[m + NUMNEIGH*(n + S)];
1050  cr[m] += Temp * (v[0] - v0[0]);
1051  cg[m] += Temp * (v[1] - v0[1]);
1052  cb[m] += Temp * (v[2] - v0[2]);
1053  }
1054  }
1055  }
1056  }
1057 
1058  /* The first three coeff values have UK_FRACBITS fractional bits */
1059  CoeffPtr[0] = v0[0] << (UK_FRACBITS - INPUT_FRACBITS);
1060  CoeffPtr[1] = v0[1] << (UK_FRACBITS - INPUT_FRACBITS);
1061  CoeffPtr[2] = v0[2] << (UK_FRACBITS - INPUT_FRACBITS);
1062 
1063  /* The other values have COEFF_FRACBITS fractional bits */
1064  for(m = 0; m < NUMNEIGH; m++)
1065  {
1066  CoeffPtr[3*(m + 1) + 0] = FLOAT_TO_FIXED(cr[m]
1068  CoeffPtr[3*(m + 1) + 1] = FLOAT_TO_FIXED(cg[m]
1070  CoeffPtr[3*(m + 1) + 2] = FLOAT_TO_FIXED(cb[m]
1072  }
1073 
1074  }
1075 
1076  for(y = 0, k = 0; y < OutputHeight; y++)
1077  {
1078  Y = YStart + (float)(y/Param.ScaleFactor);
1079  iy = (int)ceil(Y - 2*NEIGHRADIUS);
1080 
1081  /* Precompute y-factor of the window weights */
1082  for(my = 0; my < 4*NEIGHRADIUS; my++)
1083  WindowWeightY[my] = FLOAT_TO_FIXED(
1084  CubicBSpline(Y - (iy + my)), WINDOW_FRACBITS);
1085 
1086  for(x = 0; x < OutputWidth; x++, k++)
1087  {
1088  X = XStart + (float)(x/Param.ScaleFactor);
1089  ix = (int)ceil(X - 2*NEIGHRADIUS);
1090 
1091  /* Precompute x-factor of the window weights */
1092  for(mx = 0; mx < 4*NEIGHRADIUS; mx++)
1093  WindowWeightX[mx] = FLOAT_TO_FIXED(
1094  CubicBSpline(X - (ix + mx)), WINDOW_FRACBITS);
1095 
1096  DenomSum = 0;
1097  u[0] = u[1] = u[2] = 0;
1098 
1099  for(my = 0, Ypf = (int32_t)ROUND((Y - iy) * FIXED_ONE(XY_FRACBITS)); my < 4*NEIGHRADIUS;
1100  my++, Ypf -= FIXED_ONE(XY_FRACBITS))
1101  if((iy + my) >= 0 && (iy + my) < InputHeight)
1102  {
1103  for(mx = 0, Xpf = (int32_t)ROUND((X - ix) * FIXED_ONE(XY_FRACBITS)),
1104  i = ix + InputWidth*(iy + my), CoeffPtr = Coeff + i*3*(NUMNEIGH + 1);
1105  mx < 4*NEIGHRADIUS; mx++, i++, CoeffPtr += 3*(NUMNEIGH + 1), Xpf -= FIXED_ONE(XY_FRACBITS))
1106  {
1107  if((ix + mx) < 0 || (ix + mx) >= InputWidth)
1108  continue;
1109 
1110  /* WindowWeight has 2*WINDOW_FRACBITS fractional bits. */
1111  WindowWeight = (WindowWeightX[mx] * WindowWeightY[my]);
1112  /* DenomSum is computed using 2*WINDOW_FRACBITS. */
1113  DenomSum += WindowWeight;
1114  /* Now reduce to WindowWeight to WINDOW_FRACBITS. */
1115  WindowWeight = (WindowWeight + FIXED_HALF(WINDOW_FRACBITS))
1116  >> WINDOW_FRACBITS;
1117 
1118  if(!WindowWeight)
1119  continue;
1120 
1121  S = Stencil[i];
1122 
1123  uk[0] = CoeffPtr[0];
1124  uk[1] = CoeffPtr[1];
1125  uk[2] = CoeffPtr[2];
1126 
1127  for(ny = -NEIGHRADIUS, n = 3; ny <= NEIGHRADIUS; ny++)
1128  for(nx = -NEIGHRADIUS; nx <= NEIGHRADIUS; nx++, n += 3)
1129  {
1130  int32_t phit, phin;
1131 
1132  /* The tables use TRIG_FRACBITS fractional bits,
1133  and X and Y use XY_FRACBITS, so the products
1134  have (TRIG_FRACBITS + XY_FRACBITS) fractional
1135  bits. The shift reduces the result to
1136  PHITN_FRACBITS. */
1137  phit = ( CosTableTf[S]*(Xpf - nx*FIXED_ONE(XY_FRACBITS))
1138  + SinTableTf[S]*(Ypf - ny*FIXED_ONE(XY_FRACBITS))
1141  phin = ( -SinTableNf[S]*(Xpf - nx*FIXED_ONE(XY_FRACBITS))
1142  + CosTableNf[S]*(Ypf - ny*FIXED_ONE(XY_FRACBITS))
1145 
1146  /* phit and phin have PHITN_FRACBITS, so the
1147  products have 2*PHITN_FRACBITS. The result is
1148  shifted by 2*PHITN_FRACBITS to convert to
1149  quantity (by floor rounding) to an integer. */
1150  Cur = (phit*phit + phin*phin) >> (2*PHITN_FRACBITS);
1151 
1152  if(Cur >= ExpTableSize)
1153  continue;
1154 
1155  /* Compute exp(-Cur) via table look up. The result
1156  has PHI_FRACBITS fractional bits. */
1157  Weight = ExpTable[Cur];
1158 
1159  /* The Coeff values have COEFF_FRACBITS fractional
1160  bits and Weight has PHI_FRACBITS. The products
1161  are shifted so that the result has
1162  WINDOW_FRACBITS. */
1163  uk[0] += (CoeffPtr[n + 0] * Weight
1166  uk[1] += (CoeffPtr[n + 1] * Weight
1169  uk[2] += (CoeffPtr[n + 2] * Weight
1172  }
1173 
1174  /* u is computed using WINDOW_FRACBITS + UK_FRACBITS. */
1175  u[0] += WindowWeight * uk[0];
1176  u[1] += WindowWeight * uk[1];
1177  u[2] += WindowWeight * uk[2];
1178  }
1179  }
1180 
1181  if(DenomSum >= FIXED_ONE(2*WINDOW_FRACBITS) - 1)
1182  {
1183  u[0] = ROUND_FIXED(u[0], WINDOW_FRACBITS + UK_FRACBITS);
1184  u[1] = ROUND_FIXED(u[1], WINDOW_FRACBITS + UK_FRACBITS);
1185  u[2] = ROUND_FIXED(u[2], WINDOW_FRACBITS + UK_FRACBITS);
1186  }
1187  else
1188  {
1189  /* Reduce DenomSum from 2*WINDOW_FRACBITS fractional bits to
1190  (WINDOW_FRACBITS + UK_FRACBITS) fractional bits. */
1191 #if WINDOW_FRACBITS > UK_FRACBITS
1192  DenomSum = (DenomSum + FIXED_HALF(WINDOW_FRACBITS - UK_FRACBITS))
1194 #endif
1195  /* u and DenomSum both have (WINDOW_FRACBITS + UK_FRACBITS)
1196  fractional bits, so the quotient is integer. */
1197  u[0] = (u[0] + DenomSum/2) / DenomSum;
1198  u[1] = (u[1] + DenomSum/2) / DenomSum;
1199  u[2] = (u[2] + DenomSum/2) / DenomSum;
1200  }
1201 
1202  Pixel = 0xFFFFFFFF;
1203  ((uint8_t *)&Pixel)[0] = CLAMP(u[0],0,255);
1204  ((uint8_t *)&Pixel)[1] = CLAMP(u[1],0,255);
1205  ((uint8_t *)&Pixel)[2] = CLAMP(u[2],0,255);
1206  Output[k] = Pixel;
1207  }
1208  }
1209 
1210  Success = 1;
1211 Catch:
1212  Free(ExpTable);
1213  Free(Coeff);
1214  return Success;
1215 }
1216 
1217 
1219 static void AddResidual(int32_t *InputAdjusted, const int32_t *Residual,
1220  int InputWidth, int InputHeight, int PadInput)
1221 {
1222  const int PadWidth = InputWidth + 2*PadInput;
1223  const int RowEl = 3*InputWidth;
1224  const int PadRowEl = 3*PadWidth;
1225  int i, Row;
1226 
1227  Residual += 3*(PadInput + PadInput*PadWidth);
1228 
1229  for(Row = InputHeight; Row; Row--)
1230  {
1231  for(i = 0; i < RowEl; i++)
1232  InputAdjusted[i] += Residual[i];
1233 
1234  InputAdjusted += RowEl;
1235  Residual += PadRowEl;
1236  }
1237 }
1238 
1239 
1241 static void StencilStripPad(int *Stencil,
1242  int InputWidth, int InputHeight, int PadInput, int StencilMul)
1243 {
1244  const int PadWidth = InputWidth + 2*PadInput;
1245  int *Src;
1246  int i, Row;
1247 
1248  Src = Stencil + (PadInput + PadInput*PadWidth);
1249 
1250  for(Row = InputHeight; Row; Row--)
1251  {
1252  for(i = 0; i < InputWidth; i++)
1253  Stencil[i] = Src[i] / StencilMul;
1254 
1255  Stencil += InputWidth;
1256  Src += PadWidth;
1257  }
1258 }
1259 
1270 int CWInterpEx(uint32_t *Output, int OutputWidth, int OutputHeight,
1271  const uint32_t *Input, int InputWidth, int InputHeight,
1272  const int32_t *Psi, cwparams Param)
1273 {
1274  const int ScaleFactor = (int)ceil(Param.ScaleFactor);
1275  const int SupportRadius = (NEIGHRADIUS+1)*ScaleFactor - 1;
1276  const int SupportWidth = 2*SupportRadius + 1;
1277  const int SupportSize = SupportWidth*SupportWidth;
1278  const int StencilMul = NUMNEIGH*SupportSize;
1279  double *InverseA = NULL;
1280  int *Stencil = NULL;
1281  int32_t *InputFixed = NULL, *InputAdjusted = NULL, *OutputFixed = NULL, *Residual = NULL;
1282  unsigned long StartTime, StopTime;
1283  int i, PadInput, pw, ph, Success = 0;
1284  int32_t ResNorm;
1285 
1286 
1287  /* Iterative refinement is unnecessary if PSF is the Dirac delta */
1288  if(Param.PsfSigma == 0.0)
1289  Param.RefinementSteps = 0;
1290 
1291  PadInput = 4 + (ScaleFactor + 1)/2;
1292  pw = InputWidth + 2*PadInput;
1293  ph = InputHeight + 2*PadInput;
1294 
1295  if( !(InverseA = (double *)Malloc(sizeof(double)*NUMNEIGH*NUMNEIGH*NUMSTENCILS))
1296  || !(OutputFixed = (int32_t *)Malloc(sizeof(int32_t)*
1297  PIXEL_STRIDE*pw*ScaleFactor*ph*ScaleFactor))
1298  || !(InputFixed = (int32_t *)Malloc(sizeof(int32_t)*PIXEL_STRIDE*pw*ph))
1299  || !(InputAdjusted = (int32_t *)Malloc(sizeof(int32_t)*PIXEL_STRIDE*InputWidth*InputHeight))
1300  || !(Residual = (int32_t *)Malloc(sizeof(int32_t)*PIXEL_STRIDE*pw*ph))
1301  || !(Stencil = (int *)Malloc(sizeof(int)*pw*ph)) )
1302  goto Catch;
1303 
1304  if(!ComputeMatrices(InverseA, Param))
1305  goto Catch;
1306 
1307  /* Start timing */
1308  StartTime = Clock();
1309 
1310  if(Param.RefinementSteps > 0)
1311  {
1312  /* Convert 32-bit RGBA pixels to integer array */
1313  ConvertInput(InputFixed, Input, InputWidth, InputHeight, PadInput);
1314 
1315  memset(InputAdjusted, 0, sizeof(int32_t)*3*InputWidth*InputHeight);
1316  AddResidual(InputAdjusted, InputFixed, InputWidth, InputHeight, PadInput);
1317 
1318  /* Select the best-fitting contour stencils */
1319  if(!FitStencils(Stencil, InputFixed, pw, ph, StencilMul))
1320  goto Catch;
1321 
1322  memset(OutputFixed, 0, sizeof(int32_t)*
1323  3*pw*ScaleFactor*ph*ScaleFactor);
1324  memset(Residual, 0, sizeof(int32_t)*3*pw*ph);
1325 
1326  printf("\n Iteration Residual norm\n -------------------------\n");
1327 
1328  /* First interpolation pass */
1329  CWFirstPass(OutputFixed, ScaleFactor, InputFixed, pw, ph, Stencil, Psi);
1330 
1331  /* Iterative refinement */
1332  for(i = 1; i <= Param.RefinementSteps; i++)
1333  {
1334  /* Compute the residual */
1335  if((ResNorm = CWResidual(Residual, OutputFixed, InputFixed,
1336  pw, ph, Param)) < 0.0)
1337  goto Catch;
1338 
1339  printf(" %8d %15.8f\n", i, ResNorm/(255.0*256.0));
1340 
1341  AddResidual(InputAdjusted, Residual, InputWidth, InputHeight, PadInput);
1342 
1343  if(i < Param.RefinementSteps)
1344  {
1345  /* Interpolation refinement pass */
1346  CWRefinementPass(OutputFixed, ScaleFactor, Residual, pw, ph, Stencil, Psi);
1347  }
1348  }
1349 
1350  StencilStripPad(Stencil, InputWidth, InputHeight, PadInput, StencilMul);
1351  }
1352  else
1353  {
1354  /* Convert 32-bit RGBA pixels to integer array */
1355  ConvertInput(InputAdjusted, Input, InputWidth, InputHeight, 0);
1356 
1357  /* Select the best-fitting contour stencils */
1358  if(!FitStencils(Stencil, InputAdjusted, InputWidth, InputHeight, 1))
1359  goto Catch;
1360  }
1361 
1362  if(!CWArbitraryInterp(Output, OutputWidth, OutputHeight,
1363  InputAdjusted, InputWidth, InputHeight, Stencil, InverseA, Param))
1364  goto Catch;
1365 
1366  /* The final interpolation is now complete, stop timing. */
1367  StopTime = Clock();
1368 
1369  if(Param.RefinementSteps > 1)
1370  printf(" %8d (not computed)\n\n", Param.RefinementSteps + 1);
1371 
1372  /* Display the CPU time spent performing the interpolation. */
1373  printf(" CPU time: %.3f s\n\n", 0.001*(StopTime - StartTime));
1374 
1375  Success = 1;
1376 
1377 Catch: /* This label is used for error handling. If something went wrong
1378  above (which may be out of memory or a computation error), then
1379  execution jumps to this point to clean up and exit. */
1380  Free(Stencil);
1381  Free(Residual);
1382  Free(InputAdjusted);
1383  Free(InputFixed);
1384  Free(OutputFixed);
1385  Free(InverseA);
1386  return Success;
1387 }
1388 
1389 
1391 int DisplayContours(uint32_t *Output, int OutputWidth, int OutputHeight,
1392  uint32_t *Input, int InputWidth, int InputHeight, cwparams Param)
1393 {
1394  const int Pad = 2;
1395  const float LineColor[3] = {0, 0, 0};
1396  int *Stencil = 0;
1397  float dx, dy;
1398  int32_t *InputInt = NULL;
1399  uint32_t Pixel;
1400  int x, y, S, pw, ph, Success = 0;
1401 
1402 
1403  pw = InputWidth + 2*Pad;
1404  ph = InputHeight + 2*Pad;
1405 
1406  if( !(InputInt = (int32_t *)Malloc(sizeof(int32_t)*3*pw*ph))
1407  || !(Stencil = (int *)Malloc(sizeof(int)*pw*ph)) )
1408  goto Catch;
1409 
1410  /* Convert 32-bit RGBA pixels to integer array */
1411  ConvertInput(InputInt, Input, InputWidth, InputHeight, Pad);
1412 
1413  /* Select the best-fitting contour stencils */
1414  if(!FitStencils(Stencil, InputInt, pw, ph, 1))
1415  goto Catch;
1416 
1417  /* Lighten the image */
1418  for(y = 0; y < InputHeight; y++)
1419  for(x = 0; x < InputWidth; x++)
1420  {
1421  Pixel = Input[x + InputWidth*y];
1422  ((uint8_t*)&Pixel)[0] = (uint8_t)(((uint8_t*)&Pixel)[0]/2 + 128);
1423  ((uint8_t*)&Pixel)[1] = (uint8_t)(((uint8_t*)&Pixel)[1]/2 + 128);
1424  ((uint8_t*)&Pixel)[2] = (uint8_t)(((uint8_t*)&Pixel)[2]/2 + 128);
1425  Input[x + InputWidth*y] = Pixel;
1426  }
1427 
1428  /* Nearest neighbor interpolation */
1429  NearestInterp(Output, OutputWidth, OutputHeight,
1430  Input, InputWidth, InputHeight,
1431  (float)Param.ScaleFactor, Param.CenteredGrid);
1432 
1433  /* Draw contour orientation lines */
1434  for(y = 0; y < InputHeight; y++)
1435  for(x = 0; x < InputWidth; x++)
1436  {
1437  S = Stencil[(x + Pad) + pw*(y + Pad)];
1438  dx = (float)cos(StencilOrientation[S])*0.6f;
1439  dy = (float)sin(StencilOrientation[S])*0.6f;
1440  DrawLine(Output, OutputWidth, OutputHeight,
1441  (float)Param.ScaleFactor*(x - dx + 0.5f) - 0.5f,
1442  (float)Param.ScaleFactor*(y - dy + 0.5f) - 0.5f,
1443  (float)Param.ScaleFactor*(x + dx + 0.5f) - 0.5f,
1444  (float)Param.ScaleFactor*(y + dy + 0.5f) - 0.5f, LineColor);
1445  }
1446 
1447  Success = 1;
1448 Catch: /* This label is used for error handling. If something went wrong
1449  above (which may be out of memory or a computation error), then
1450  execution jumps to this point to clean up and exit. */
1451  Free(Stencil);
1452  Free(InputInt);
1453  return Success;
1454 }
1455 
1456