22 #include "mosaic_util.h"
24 #include "create_xgrid.h"
50 void cubic_spline_sp(
int size1,
int size2,
const double *grid1,
const double *grid2,
const double *data1,
53 double *delta=NULL, *d=NULL, *dh=NULL, *b=NULL, *c = NULL;
55 int i, k, n, klo, khi, kmax;
57 for(i=1; i<size1; i++) {
58 if( grid1[i] <= grid1[i-1] ) error_handler(
"cubic_spline_sp: grid1 is not monotonic increasing");
61 for(i=0; i<size2; i++) {
62 if( grid2[i] < grid1[0] || grid2[i] > grid1[size1-1]) error_handler(
"cubic_spline_sp: grid2 lies outside grid1");
65 if(size1 < 2) error_handler(
"cubic_spline_sp: the size of input grid should be at least 2");
67 p = (data1[1]-data1[0])/(grid1[1]-grid1[0]);
68 for(i=0; i< size2; i++) data2[i] = p*(grid2[i] - grid1[0]) + data1[0];
71 delta = (
double *)malloc((size1-1)*
sizeof(double));
72 dh = (
double *)malloc((size1-1)*
sizeof(double));
73 d = (
double *)malloc(size1*
sizeof(
double));
74 for(k=0;k<size1-1;k++) {
75 dh[k] = grid1[k+1]-grid1[k];
76 delta[k] = (data1[k+1]-data1[k])/(dh[k]);
81 for(k=1;k<size1-1;k++) {
82 if( delta[k]*delta[k-1] > 0.0 ) {
83 w1 = 2.0*dh[k] + dh[k-1];
84 w2 = dh[k] + 2.0*dh[k-1];
85 d[k] = (w1+w2)/(w1/delta[k-1]+w2/delta[k]);
95 d[0] = ((2.0*dh[0] + dh[1])*delta[0] - dh[0]*delta[1])/(dh[0]+dh[1]);
97 if ( d[0]*delta[0] < 0.0 ) {
101 if ( delta[0]*delta[1] < 0.0 && fabs(d[0]) > fabs(3.0*delta[0])) {
106 d[kmax] = ((2.0*dh[kmax-1] + dh[kmax-2])*delta[kmax-1] - dh[kmax-1]*delta[kmax-2])/(dh[kmax-1]+dh[kmax-2]);
107 if ( d[kmax]*delta[kmax-1] < 0.0 ) {
111 if ( delta[kmax-1]*delta[kmax-2] < 0.0 && fabs(d[kmax]) > fabs(3.0*delta[kmax-1])) {
112 d[kmax]=3.0*delta[kmax-1];
117 b = (
double *)malloc((size1-1)*
sizeof(double));
118 c = (
double *)malloc((size1-1)*
sizeof(double));
119 for (k=0; k<size1-1; k++) {
120 c[k] = (3.0*delta[k]-2.0*d[k]-d[k+1])/dh[k];
121 b[k] = (d[k]-2.0*delta[k]+d[k+1])/(dh[k]*dh[k]);
124 for(k=0; k<size2; k++) {
125 n = nearest_index(grid2[k],grid1, size1);
126 if (grid1[n] < grid2[k]) {
138 s = grid2[k] - grid1[klo];
139 data2[k] = data1[klo]+s*(d[klo]+s*(c[klo]+s*b[klo]));
175 void cubic_spline(
int size1,
int size2,
const double *grid1,
const double *grid2,
const double *data1,
176 double *data2,
double yp1,
double ypn )
178 double *y2=NULL, *u=NULL;
179 double p, sig, qn, un, h, a, b;
180 int i, k, n, klo, khi;
182 for(i=1; i<size1; i++) {
183 if( grid1[i] <= grid1[i-1] ) error_handler(
"cubic_spline: grid1 is not monotonic increasing");
186 for(i=0; i<size2; i++) {
187 if( grid2[i] < grid1[0] || grid2[i] > grid1[size1-1]) error_handler(
"cubic_spline: grid2 lies outside grid1");
190 if(size1 < 2) error_handler(
"cubic_spline: the size of input grid should be at least 2");
192 p = (data1[1]-data1[0])/(grid1[1]-grid1[0]);
193 for(i=0; i< size2; i++) data2[i] = p*(grid2[i] - grid1[0]) + data1[0];
196 y2 = (
double *)malloc(size1*
sizeof(
double));
197 u = (
double *)malloc(size1*
sizeof(
double));
204 u[0]=(3./(grid1[1]-grid1[0]))*((data1[1]-data1[0])/(grid1[1]-grid1[0])-yp1);
207 for(i=1; i<size1-1; i++) {
208 sig=(grid1[i]-grid1[i-1])/(grid1[i+1]-grid1[i-1]);
211 u[i]=(6.*((data1[i+1]-data1[i])/(grid1[i+1]-grid1[i])-(data1[i]-data1[i-1])
212 /(grid1[i]-grid1[i-1]))/(grid1[i+1]-grid1[i-1])-sig*u[i-1])/p;
222 un=(3./(grid1[size1-1]-grid1[size1-2]))*(ypn-(data1[size1-1]-data1[size1-2])/(grid1[size1-1]-grid1[size1-2]));
225 y2[size1-1]=(un-qn*u[size1-2])/(qn*y2[size1-2]+1.);
227 for(k=size1-2; k>=0; k--) y2[k] = y2[k]*y2[k+1]+u[k];
230 for(k=0; k<size2; k++) {
231 n = nearest_index(grid2[k],grid1, size1);
232 if (grid1[n] < grid2[k]) {
244 h = grid1[khi]-grid1[klo];
245 a = (grid1[khi] - grid2[k])/h;
246 b = (grid2[k] - grid1[klo])/h;
247 data2[k] = a*data1[klo] + b*data1[khi]+ ((pow(a,3.0)-a)*y2[klo] + (pow(b,3.0)-b)*y2[khi])*(pow(h,2.0))/6;
261 void conserve_interp(
int nx_src,
int ny_src,
int nx_dst,
int ny_dst,
const double *x_src,
262 const double *y_src,
const double *x_dst,
const double *y_dst,
263 const double *mask_src,
const double *data_src,
double *data_dst )
266 int *xgrid_i1, *xgrid_j1, *xgrid_i2, *xgrid_j2;
267 double *xgrid_area, *dst_area, *area_frac;
270 xgrid_i1 = (
int *)malloc(MAXXGRID*
sizeof(
int));
271 xgrid_j1 = (
int *)malloc(MAXXGRID*
sizeof(
int));
272 xgrid_i2 = (
int *)malloc(MAXXGRID*
sizeof(
int));
273 xgrid_j2 = (
int *)malloc(MAXXGRID*
sizeof(
int));
274 xgrid_area = (
double *)malloc(MAXXGRID*
sizeof(
double));
275 dst_area = (
double *)malloc(nx_dst*ny_dst*
sizeof(
double));
276 nxgrid = create_xgrid_2dx2d_order1(&nx_src, &ny_src, &nx_dst, &ny_dst, x_src, y_src, x_dst, y_dst, mask_src,
277 xgrid_i1, xgrid_j1, xgrid_i2, xgrid_j2, xgrid_area );
282 for(n=0; n<nx_dst*ny_dst; n++) dst_area[n] = 0;
283 for(n=0; n<nxgrid; n++) dst_area[xgrid_j2[n]*nx_dst+xgrid_i2[n]] += xgrid_area[n];
285 area_frac = (
double *)malloc(nxgrid*
sizeof(
double));
286 for(n=0; n<nxgrid; n++) area_frac[n] = xgrid_area[n]/dst_area[xgrid_j2[n]*nx_dst+xgrid_i2[n]];
288 for(n=0; n<nx_dst*ny_dst; n++) {
291 for(n=0; n<nxgrid; n++) {
292 data_dst[xgrid_j2[n]*nx_dst+xgrid_i2[n]] += data_src[xgrid_j1[n]*nx_src+xgrid_i1[n]]*area_frac[n];
311 void conserve_interp_great_circle(
int nx_src,
int ny_src,
int nx_dst,
int ny_dst,
const double *x_src,
312 const double *y_src,
const double *x_dst,
const double *y_dst,
313 const double *mask_src,
const double *data_src,
double *data_dst )
316 int *xgrid_i1, *xgrid_j1, *xgrid_i2, *xgrid_j2;
317 double *xgrid_area, *dst_area, *area_frac, *xgrid_di, *xgrid_dj;
320 xgrid_i1 = (
int *)malloc(MAXXGRID*
sizeof(
int));
321 xgrid_j1 = (
int *)malloc(MAXXGRID*
sizeof(
int));
322 xgrid_i2 = (
int *)malloc(MAXXGRID*
sizeof(
int));
323 xgrid_j2 = (
int *)malloc(MAXXGRID*
sizeof(
int));
324 xgrid_area = (
double *)malloc(MAXXGRID*
sizeof(
double));
325 xgrid_di = (
double *)malloc(MAXXGRID*
sizeof(
double));
326 xgrid_dj = (
double *)malloc(MAXXGRID*
sizeof(
double));
327 dst_area = (
double *)malloc(nx_dst*ny_dst*
sizeof(
double));
328 nxgrid = create_xgrid_great_circle(&nx_src, &ny_src, &nx_dst, &ny_dst, x_src, y_src, x_dst, y_dst, mask_src,
329 xgrid_i1, xgrid_j1, xgrid_i2, xgrid_j2, xgrid_area, xgrid_di, xgrid_dj );
334 for(n=0; n<nx_dst*ny_dst; n++) dst_area[n] = 0;
335 for(n=0; n<nxgrid; n++) dst_area[xgrid_j2[n]*nx_dst+xgrid_i2[n]] += xgrid_area[n];
337 area_frac = (
double *)malloc(nxgrid*
sizeof(
double));
338 for(n=0; n<nxgrid; n++) area_frac[n] = xgrid_area[n]/dst_area[xgrid_j2[n]*nx_dst+xgrid_i2[n]];
340 for(n=0; n<nx_dst*ny_dst; n++) {
343 for(n=0; n<nxgrid; n++) {
344 data_dst[xgrid_j2[n]*nx_dst+xgrid_i2[n]] += data_src[xgrid_j1[n]*nx_src+xgrid_i1[n]]*area_frac[n];
359 void linear_vertical_interp(
int nx,
int ny,
int nk1,
int nk2,
const double *grid1,
const double *grid2,
360 double *data1,
double *data2)
365 for(k=1; k<nk1; k++) {
366 if(grid1[k] <= grid1[k-1]) error_handler(
"interp.c: grid1 not monotonic");
368 for(k=1; k<nk2; k++) {
369 if(grid2[k] <= grid2[k-1]) error_handler(
"interp.c: grid2 not monotonic");
372 if (grid1[0] > grid2[0] ) error_handler(
"interp.c: grid2 lies outside grid1");
373 if (grid1[nk1-1] < grid2[nk2-1] ) error_handler(
"interp.c: grid2 lies outside grid1");
375 for(k=0; k<nk2; k++) {
376 n = nearest_index(grid2[k],grid1,nk1);
377 if (grid1[n] < grid2[k]) {
378 w = (grid2[k]-grid1[n])/(grid1[n+1]-grid1[n]);
379 for(l=0; l<nx*ny; l++) {
380 data2[k*nx*ny+l] = (1.-w)*data1[n*nx*ny+l] + w*data1[(n+1)*nx*ny+l];
385 for(l=0;l<nx*ny;l++) data2[k*nx*ny+l] = data1[n*nx*ny+l];
387 w = (grid2[k]-grid1[n-1])/(grid1[n]-grid1[n-1]);
388 for(l=0; l<nx*ny; l++) {
389 data2[k*nx*ny+l] = (1.-w)*data1[(n-1)*nx*ny+l] + w*data1[n*nx*ny+l];