//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_set_ri(HADT d, HCDT re, HCDT im) { long nx=d->GetNx(),ny=d->GetNy(),nz=d->GetNz(); d->Create(nx,ny,nz); #pragma omp parallel for for(long i=0;i<nx*ny*nz;i++) d->a[i] = dual(re->vthr(i),im->vthr(i)); }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_set_ap(HADT d, HCDT a, HCDT p) { long nx=d->GetNx(),ny=d->GetNy(),nz=d->GetNz(); d->Create(nx,ny,nz); #pragma omp parallel for for(long i=0;i<nx*ny*nz;i++) { register mreal aa=a->vthr(i), pp=p->vthr(i); d->a[i] = dual(aa*cos(pp), aa*sin(pp)); } }
//----------------------------------------------------------------------------- int MGL_EXPORT mgl_datac_read_all(HADT dat, const char *templ, int as_slice) { #ifndef WIN32 mglDataC d; glob_t res; size_t i; dual *b; long kx,ky,kz; glob (templ, GLOB_TILDE, NULL, &res); //read first file for(i=0;i<res.gl_pathc;i++) if(mgl_datac_read(&d,res.gl_pathv[i])) break; if(i>=res.gl_pathc) { globfree (&res); return false; } kx = d.nx; ky = d.ny; kz = d.nz; b = (dual *)malloc(kx*ky*kz*sizeof(dual)); memcpy(b,d.a,kx*ky*kz*sizeof(dual)); for(;i<res.gl_pathc;i++) { if(mgl_datac_read(&d,res.gl_pathv[i])) if(!mgl_add_file(kx,ky,kz,b,&d,as_slice)) { globfree (&res); free(b); return false; } } dat->Set(b,kx,ky,kz); globfree (&res); free(b); return true; #else return false; #endif }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_squeeze(HADT d, long rx,long ry,long rz,long smooth) { long kx,ky,kz, nx=d->nx, ny=d->ny, nz=d->nz; dual *b; // simple checking if(rx>=nx) rx=nx-1; if(rx<1) rx=1; if(ry>=ny) ry=ny-1; if(ry<1) ry=1; if(rz>=nz) rz=nz-1; if(rz<1) rz=1; // new sizes kx = 1+(nx-1)/rx; ky = 1+(ny-1)/ry; kz = 1+(nz-1)/rz; b = new dual[kx*ky*kz]; if(!smooth) #pragma omp parallel for collapse(3) for(long k=0;k<kz;k++) for(long j=0;j<ky;j++) for(long i=0;i<kx;i++) b[i+kx*(j+ky*k)] = d->a[i*rx+nx*(j*ry+ny*rz*k)]; else #pragma omp parallel for collapse(3) for(long k=0;k<kz;k++) for(long j=0;j<ky;j++) for(long i=0;i<kx;i++) { long dx,dy,dz,i1,j1,k1; dx = (i+1)*rx<=nx ? rx : nx-i*rx; dy = (j+1)*ry<=ny ? ry : ny-j*ry; dz = (k+1)*rz<=nz ? rz : nz-k*rz; dual s = 0; for(k1=k*rz;k1<k*rz+dz;k1++) for(j1=j*ry;j1<j*ry+dz;j1++) for(i1=i*rx;i1<i*rx+dx;i1++) s += d->a[i1+nx*(j1+ny*k1)]; b[i+kx*(j+ky*k)] = s/mreal(dx*dy*dz); } if(!d->link) delete [](d->a); d->a=b; d->nx = kx; d->ny = ky; d->nz = kz; d->NewId(); d->link=false; }
//----------------------------------------------------------------------------- int MGL_EXPORT mgl_datac_read_range(HADT dat, const char *templ, double from, double to, double step, int as_slice) { mglDataC d; double t = from; dual *b; long kx,ky,kz,n=strlen(templ)+20; char *fname = new char[n]; //read first file do{ snprintf(fname,n,templ,t); t+= step; } while(!mgl_datac_read(&d,fname) && t<=to); if(t>to) { delete []fname; return false; } kx = d.nx; ky = d.ny; kz = d.nz; b = (dual *)malloc(kx*ky*kz*sizeof(dual)); memcpy(b,d.a,kx*ky*kz*sizeof(dual)); // read other files for(;t<=to;t+=step) { snprintf(fname,n,templ,t); if(mgl_datac_read(&d,fname)) if(!mgl_add_file(kx,ky,kz,b,&d,as_slice)) { delete []fname; free(b); return false; } } dat->Set(b,kx,ky,kz); delete []fname; free(b); return true; }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_link(HADT d, dual *A, long mx,long my,long mz) { if(!A) return; if(!d->link && d->a) delete [](d->a); d->nx = mx>0 ? mx:1; d->ny = my>0 ? my:1; d->nz = mz>0 ? mz:1; d->link=true; d->a=A; d->NewId(); }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_rearrange(HADT d, long mx,long my,long mz) { if(mx<1) return; // wrong mx if(my<1) { my = d->nx*d->ny*d->nz/mx; mz = 1; } else if(mz<1) mz = (d->nx*d->ny*d->nz)/(mx*my); long m = mx*my*mz; if(m==0 || m>d->nx*d->ny*d->nz) return; // too high desired dimensions d->nx = mx; d->ny = my; d->nz = mz; d->NewId(); }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_fill_eq(HMGL gr, HADT d, const char *eq, HCDT vdat, HCDT wdat, const char *opt) { if(vdat && vdat->GetNN()!=d->GetNN()) return; // incompatible dimensions if(wdat && wdat->GetNN()!=d->GetNN()) return; gr->SaveState(opt); std::wstring s = d->Name(); d->Name(L"u"); mglDataV x(d->nx,d->ny,d->nz, gr->Min.x,gr->Max.x,'x'); x.Name(L"x"); mglDataV y(d->nx,d->ny,d->nz, gr->Min.y,gr->Max.y,'y'); y.Name(L"y"); mglDataV z(d->nx,d->ny,d->nz, gr->Min.z,gr->Max.z,'z'); z.Name(L"z"); mglDataV i(d->nx,d->ny,d->nz, 0,d->nx-1,'x'); i.Name(L"i"); mglDataV j(d->nx,d->ny,d->nz, 0,d->ny-1,'y'); j.Name(L"j"); mglDataV k(d->nx,d->ny,d->nz, 0,d->nz-1,'z'); k.Name(L"k"); mglDataV r(d->nx,d->ny,d->nz); r.Name(L"#$mgl"); mglData v(vdat), w(wdat); v.Name(L"v"); w.Name(L"w"); std::vector<mglDataA*> list; list.push_back(&x); list.push_back(&y); list.push_back(&z); list.push_back(&r); list.push_back(d); list.push_back(&v); list.push_back(&w); list.push_back(&i); list.push_back(&j); list.push_back(&k); d->Move(mglFormulaCalcC(eq,list)); d->Name(s.c_str()); gr->LoadState(); }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_transpose(HADT d, const char *dim) { long nx=d->nx, ny=d->ny, nz=d->nz, n; dual *b=new dual[nx*ny*nz], *a=d->a; if(!strcmp(dim,"xyz")) memcpy(b,a,nx*ny*nz*sizeof(dual)); else if(!strcmp(dim,"xzy") || !strcmp(dim,"zy")) { #pragma omp parallel for collapse(3) for(long j=0;j<ny;j++) for(long k=0;k<nz;k++) for(long i=0;i<nx;i++) b[i+nx*(k+nz*j)] = a[i+nx*(j+ny*k)]; n=nz; nz=ny; ny=n; } else if(!strcmp(dim,"yxz") || !strcmp(dim,"yx")) { #pragma omp parallel for collapse(3) for(long k=0;k<nz;k++) for(long i=0;i<nx;i++) for(long j=0;j<ny;j++) b[j+ny*(i+nx*k)] = a[i+nx*(j+ny*k)]; n=nx; nx=ny; ny=n; } else if(!strcmp(dim,"yzx")) { #pragma omp parallel for collapse(3) for(long k=0;k<nz;k++) for(long i=0;i<nx;i++) for(long j=0;j<ny;j++) b[j+ny*(k+nz*i)] = a[i+nx*(j+ny*k)]; n=nx; nx=ny; ny=nz; nz=n; } else if(!strcmp(dim,"zxy")) { #pragma omp parallel for collapse(3) for(long i=0;i<nx;i++) for(long j=0;j<ny;j++) for(long k=0;k<nz;k++) b[k+nz*(i+nx*j)] = a[i+nx*(j+ny*k)]; n=nx; nx=nz; nz=ny; ny=n; } else if(!strcmp(dim,"zyx") || !strcmp(dim,"zx")) { #pragma omp parallel for collapse(3) for(long i=0;i<nx;i++) for(long j=0;j<ny;j++) for(long k=0;k<nz;k++) b[k+nz*(j+ny*i)] = a[i+nx*(j+ny*k)]; n=nz; nz=nx; nx=n; } memcpy(a,b,nx*ny*nz*sizeof(dual)); delete []b; n=d->nx; d->nx=nx; d->ny=ny; d->nz=nz; if(nx!=n) d->NewId(); }
//----------------------------------------------------------------------------- void MGL_EXPORT mgl_datac_extend(HADT d, long n1, long n2) { long nx=d->nx, ny=d->ny, nz=d->nz; if(nz>2 || n1==0) return; long mx, my, mz; dual *b=0; if(n1>0) // extend to higher dimension(s) { n2 = n2>0 ? n2:1; mx = nx; my = ny>1?ny:n1; mz = ny>1 ? n1 : n2; b = new dual[mx*my*mz]; if(ny>1) #pragma omp parallel for for(long i=0;i<n1;i++) memcpy(b+i*nx*ny, d->a, nx*ny*sizeof(dual)); else #pragma omp parallel for for(long i=0;i<n1*n2;i++) memcpy(b+i*nx, d->a, nx*sizeof(dual)); } else { mx = -n1; my = n2<0 ? -n2 : nx; mz = n2<0 ? nx : ny; if(n2>0 && ny==1) mz = n2; b = new dual[mx*my*mz]; dual v; if(n2<0) #pragma omp parallel for collapse(2) for(long j=0;j<nx;j++) for(long i=0;i<mx*my;i++) b[i+mx*my*j] = d->a[j]; else #pragma omp parallel for collapse(2) for(long j=0;j<nx*ny;j++) for(long i=0;i<mx;i++) b[i+mx*j] = d->a[j]; if(n2>0 && ny==1) #pragma omp parallel for for(long i=0;i<n2;i++) memcpy(b+i*mx*my, d->a, mx*my*sizeof(dual)); } if(!d->link) delete [](d->a); d->a=b; d->nx=mx; d->ny=my; d->nz=mz; d->NewId(); d->link=false; }
//----------------------------------------------------------------------------- // convert substrings to arguments void mglParser::FillArg(mglGraph *gr, int k, std::wstring *arg, mglArg *a) { register long n; for(n=1;n<k;n++) { mglDataA *v; mglNum *f; a[n-1].type = -1; if(arg[n][0]=='|') a[n-1].type = -1; else if(arg[n][0]=='\'') // this is string (simplest case) { a[n-1].type = 1; std::wstring &w=arg[n],f; wchar_t buf[32]; long i,i1,ll=w.length(); for(i=1;i<ll;i++) { if(w[i]=='\'') { if(i==ll-1) continue; i++; i1 = i; if(w[i1]==',') i1++; if(w[i1]==0) continue; for(;i<ll && w[i]!='\'';i++); if(i>i1) { if(w[i1]=='!') { HADT d = mglFormulaCalcC(w.substr(i1+1,i-i1-(w[i]=='\''?1:0)), this, DataList); mreal di = imag(d->a[0]), dr = real(d->a[0]); if(di>0) mglprintf(buf,32,L"%g+%gi",dr,di); else if(di<0) mglprintf(buf,32,L"%g-%gi",dr,-di); // TODO use \u2212 ??? else mglprintf(buf,32,L"%g",dr); a[n-1].w += buf; delete d; } else { HMDT d = mglFormulaCalc(w.substr(i1,i-i1-(w[i]=='\''?1:0)), this, DataList); mglprintf(buf,32,L"%g",d->a[0]); a[n-1].w += buf; delete d; } } } else a[n-1].w += w[i]; } } else if(arg[n][0]=='{') { // this is temp data mglData *u=new mglData; std::wstring s = arg[n].substr(1,arg[n].length()-2); a[n-1].w = u->s = L"/*"+s+L"*/"; a[n-1].type = 0; ParseDat(gr, s, *u); a[n-1].d = u; u->temp=true; DataList.push_back(u); } else if((v = FindVar(arg[n].c_str()))!=0) // try to find normal variables (for data creation) { a[n-1].type=0; a[n-1].d=v; a[n-1].w=v->s; } else if((f = FindNum(arg[n].c_str()))!=0) // try to find normal number (for data creation) { a[n-1].type=2; a[n-1].d=0; a[n-1].v=f->d; a[n-1].c=f->c; a[n-1].w = f->s; } else if(arg[n][0]=='!') // complex array is asked { // parse all numbers and formulas by unified way HADT d = mglFormulaCalcC(arg[n].substr(1), this, DataList); if(d->GetNN()==1) { if(CheckForName(arg[n].substr(1))) { a[n-1].type = 2; a[n-1].v = d->v(0); a[n-1].c = d->a[0]; } else { a[n-1].type = 0; a[n-1].d = AddVar(arg[n].c_str()); } delete d; } else { a[n-1].w = L"/*"+arg[n]+L"*/"; d->temp=true; DataList.push_back(d); a[n-1].type = 0; a[n-1].d = d; } } else { // parse all numbers and formulas by unified way HMDT d = mglFormulaCalc(arg[n], this, DataList); if(d->GetNN()==1) { if(CheckForName(arg[n])) { a[n-1].type = 2; a[n-1].c = a[n-1].v = d->v(0); } else { a[n-1].type = 0; a[n-1].d = AddVar(arg[n].c_str()); } delete d; } else { a[n-1].w = L"/*"+arg[n]+L"*/"; d->temp=true; DataList.push_back(d); a[n-1].type = 0; a[n-1].d = d; } } } }