static void test_cov(){/*not good */ rand_t rstat; int seed=4; double r0=0.2; double dx=1./64; long N=1+1024; long nx=N; long ny=N; long nframe=1; seed_rand(&rstat, seed); map_t *atm=mapnew(nx, ny, dx,dx, NULL); cmat *atmhat=cnew((N+1)*3,(N+1)*3); dmat *atmhattot=dnew((N+1)*3,(N+1)*3); //cfft2plan(atmhat,-1); //cfft2plan(atmhat, 1); dset((dmat*)atm,1); cembedd(atmhat, (dmat*)atm, 0); cfft2(atmhat, -1); cabs22d(&atmhattot, 1, atmhat, 1); ccpd(&atmhat, atmhattot); cfft2i(atmhat, 1); cfftshift(atmhat); dmat *denom=dnew((N+1)*3,(N+1)*3); dmat *cov=dnew((N+1)*3,(N+1)*3); creal2d(&denom, 0, atmhat, 1); writebin(denom, "denom.bin"); dzero(atmhattot); for(long i=0; i<nframe; i++){ info("%ld of %ld\n", i, nframe); for(long j=0; j<nx*ny; j++){ atm->p[j]=randn(&rstat); } fractal_do((dmat*)atm, dx, r0,L0,ninit); /*mapwrite(atm, "atm_%ld.bin", i); */ cembedd(atmhat, (dmat*)atm, 0); cfft2(atmhat, -1); cabs22d(&atmhattot, 1, atmhat, 1); if(i==0 || (i+1)%10==0){ dscale(atmhattot, 1./(i+1)); ccpd(&atmhat, atmhattot); writebin(atmhattot, "atm_psf_%ld.bin",i+1); cfft2i(atmhat, 1); cfftshift(atmhat); creal2d(&cov, 0, atmhat, 1); for(long k=0; k<cov->nx*cov->ny; k++){ cov->p[k]/=denom->p[k]; } writebin(cov, "atm_cov_%ld.bin",i+1); } } }
/** Closed loop simulation main loop. It calls init_simu() to initialize the simulation struct. Then calls genatm() to generate atmospheric turbulence screens. Then for every time step, it calls perfevl() to evaluate performance, wfsgrad() to do wfs measurement, reconstruct() to do tomography and DM fit, filter() to do DM command filtering. In MOAO mode, it call calls moao_recon() for MOAO DM fitting. \callgraph */ void maos_sim(){ const PARMS_T *parms=global->parms; POWFS_T *powfs=global->powfs; RECON_T *recon=global->recon; APER_T *aper=global->aper; int simend=parms->sim.end; int simstart=parms->sim.start; if(parms->sim.skysim){ save_skyc(powfs,recon,parms); } if(parms->evl.psfmean || parms->evl.psfhist){ /*compute diffraction limited PSF. Save to output directory.*/ dmat *iopdevl=dnew(aper->locs->nloc,1); ccell *psf2s=0; locfft_psf(&psf2s, aper->embed, iopdevl, parms->evl.psfsize, 0); const int nwvl=parms->evl.nwvl; dcell *evlpsfdl=dcellnew(nwvl,1); for(int iwvl=0; iwvl<nwvl; iwvl++){ cabs22d(&evlpsfdl->p[iwvl], 1, psf2s->p[iwvl], 1); evlpsfdl->p[iwvl]->header=evl_header(parms, aper, -1, iwvl); } ccellfree(psf2s); writebin(evlpsfdl, "evlpsfdl.fits"); dcellfree(evlpsfdl); dfree(iopdevl); } info2("PARALLEL=%d\n", PARALLEL); if(simstart>=simend) return; double restot=0; long rescount=0; for(int iseed=0; iseed<parms->sim.nseed; iseed++){ SIM_T *simu=0; while(!(simu=maos_iseed(iseed))){ iseed++; } #ifdef HAVE_NUMA_H numa_set_localalloc(); #endif for(int isim=simstart; isim<simend; isim++){ maos_isim(isim); if(parms->sim.pause){ mypause(); } }/*isim */ { /*Compute average performance*/ int isim0; if(parms->sim.closeloop){ if(parms->sim.end>100){ isim0=MAX(50,parms->sim.end/10); }else{ isim0=MIN(20, parms->sim.end/2); } }else{ isim0=0; } double sum=0; for(int i=isim0; i<parms->sim.end; i++){ sum+=simu->cle->p[i*parms->evl.nmod]; } restot+=sum/(parms->sim.end-isim0); rescount++; } free_simu(simu); global->simu=0; }/*seed */ printf("%g\n", sqrt(restot/rescount)*1e9); }
static void test_psd(){ rand_t rstat; int seed=4; double r0=0.2; double dx=1./64; long N=1024; long nx=N; long ny=N; long ratio=1; long xskip=nx*(ratio-1)/2; long yskip=ny*(ratio-1)/2; long nframe=512; seed_rand(&rstat, seed); if(1){ map_t *atm=mapnew(nx+1, ny+1, dx,dx, NULL); cmat *hat=cnew(nx*ratio, ny*ratio); //cfft2plan(hat, -1); dmat *hattot=dnew(nx*ratio, ny*ratio); for(long i=0; i<nframe; i++){ info2("%ld of %ld\n", i, nframe); for(long j=0; j<(nx+1)*(ny+1); j++){ atm->p[j]=randn(&rstat); } fractal_do((dmat*)atm, dx, r0,L0,ninit); czero(hat); for(long iy=0; iy<ny; iy++){ for(long ix=0; ix<nx; ix++){ IND(hat,ix+xskip,iy+yskip)=IND(atm,ix,iy); } } cfftshift(hat); cfft2i(hat, -1); cabs22d(&hattot, 1, hat, 1); } dscale(hattot, 1./nframe); dfftshift(hattot); writebin(hattot, "PSD_fractal"); } { dmat *spect=turbpsd(nx, ny, dx, r0, 100, 0, 0.5); writebin(spect, "spect"); cmat *hat=cnew(nx*ratio, ny*ratio); //cfft2plan(hat, -1); dmat *hattot=dnew(nx*ratio, ny*ratio); cmat *atm=cnew(nx, ny); //cfft2plan(atm, -1); dmat *atmr=dnew(atm->nx, atm->ny); dmat *atmi=dnew(atm->nx, atm->ny); cmat* phat=hat; dmat* patmr=atmr; dmat* patmi=atmi; for(long ii=0; ii<nframe; ii+=2){ info2("%ld of %ld\n", ii, nframe); for(long i=0; i<atm->nx*atm->ny; i++){ atm->p[i]=COMPLEX(randn(&rstat), randn(&rstat))*spect->p[i]; } cfft2(atm, -1); for(long i=0; i<atm->nx*atm->ny; i++){ atmr->p[i]=creal(atm->p[i]); atmi->p[i]=cimag(atm->p[i]); } czero(hat); for(long iy=0; iy<ny; iy++){ for(long ix=0; ix<nx; ix++){ IND(phat,ix+xskip,iy+yskip)=IND(patmr,ix,iy); } } cfftshift(hat); cfft2i(hat, -1); cabs22d(&hattot, 1, hat, 1); czero(hat); for(long iy=0; iy<ny; iy++){ for(long ix=0; ix<nx; ix++){ IND(phat,ix+xskip,iy+yskip)=IND(patmi,ix,iy); } } cfftshift(hat); cfft2i(hat, -1); cabs22d(&hattot, 1, hat, 1); } dscale(hattot, 1./nframe); dfftshift(hattot); writebin(hattot, "PSD_fft"); } }
/** Time domain physical simulation. noisy: - 0: no noise at all; - 1: poisson and read out noise. - 2: only poisson noise. */ dmat *skysim_sim(dmat **mresout, const dmat *mideal, const dmat *mideal_oa, double ngsol, ASTER_S *aster, const POWFS_S *powfs, const PARMS_S *parms, int idtratc, int noisy, int phystart){ int dtratc=0; if(!parms->skyc.multirate){ dtratc=parms->skyc.dtrats->p[idtratc]; } int hasphy; if(phystart>-1 && phystart<aster->nstep){ hasphy=1; }else{ hasphy=0; } const int nmod=mideal->nx; dmat *res=dnew(6,1);/*Results. 1-2: NGS and TT modes., 3-4:On axis NGS and TT modes, 4-6: On axis NGS and TT wihtout considering un-orthogonality.*/ dmat *mreal=NULL;/*modal correction at this step. */ dmat *merr=dnew(nmod,1);/*modal error */ dcell *merrm=dcellnew(1,1);dcell *pmerrm=NULL; const int nstep=aster->nstep?aster->nstep:parms->maos.nstep; dmat *mres=dnew(nmod,nstep); dmat* rnefs=parms->skyc.rnefs; dcell *zgradc=dcellnew3(aster->nwfs, 1, aster->ngs, 0); dcell *gradout=dcellnew3(aster->nwfs, 1, aster->ngs, 0); dmat *gradsave=0; if(parms->skyc.dbg){ gradsave=dnew(aster->tsa*2,nstep); } SERVO_T *st2t=0; kalman_t *kalman=0; dcell *mpsol=0; dmat *pgm=0; dmat *dtrats=0; int multirate=parms->skyc.multirate; if(multirate){ kalman=aster->kalman[0]; dtrats=aster->dtrats; }else{ if(parms->skyc.servo>0){ const double dtngs=parms->maos.dt*dtratc; st2t=servo_new(merrm, NULL, 0, dtngs, aster->gain->p[idtratc]); pgm=aster->pgm->p[idtratc]; }else{ kalman=aster->kalman[idtratc]; } } if(kalman){ kalman_init(kalman); mpsol=dcellnew(aster->nwfs, 1); //for psol grad. } const long nwvl=parms->maos.nwvl; dcell **psf=0, **mtche=0, **ints=0; ccell *wvf=0,*wvfc=0, *otf=0; if(hasphy){ psf=mycalloc(aster->nwfs,dcell*); wvf=ccellnew(aster->nwfs,1); wvfc=ccellnew(aster->nwfs,1); mtche=mycalloc(aster->nwfs,dcell*); ints=mycalloc(aster->nwfs,dcell*); otf=ccellnew(aster->nwfs,1); for(long iwfs=0; iwfs<aster->nwfs; iwfs++){ const int ipowfs=aster->wfs[iwfs].ipowfs; const long ncomp=parms->maos.ncomp[ipowfs]; const long nsa=parms->maos.nsa[ipowfs]; wvf->p[iwfs]=cnew(ncomp,ncomp); wvfc->p[iwfs]=NULL; psf[iwfs]=dcellnew(nsa,nwvl); //cfft2plan(wvf->p[iwfs], -1); if(parms->skyc.multirate){ mtche[iwfs]=aster->wfs[iwfs].pistat->mtche[(int)aster->idtrats->p[iwfs]]; }else{ mtche[iwfs]=aster->wfs[iwfs].pistat->mtche[idtratc]; } otf->p[iwfs]=cnew(ncomp,ncomp); //cfft2plan(otf->p[iwfs],-1); //cfft2plan(otf->p[iwfs],1); ints[iwfs]=dcellnew(nsa,1); int pixpsa=parms->skyc.pixpsa[ipowfs]; for(long isa=0; isa<nsa; isa++){ ints[iwfs]->p[isa]=dnew(pixpsa,pixpsa); } } } for(int irep=0; irep<parms->skyc.navg; irep++){ if(kalman){ kalman_init(kalman); }else{ servo_reset(st2t); } dcellzero(zgradc); dcellzero(gradout); if(ints){ for(int iwfs=0; iwfs<aster->nwfs; iwfs++){ dcellzero(ints[iwfs]); } } for(int istep=0; istep<nstep; istep++){ memcpy(merr->p, PCOL(mideal,istep), nmod*sizeof(double)); dadd(&merr, 1, mreal, -1);/*form NGS mode error; */ memcpy(PCOL(mres,istep),merr->p,sizeof(double)*nmod); if(mpsol){//collect averaged modes for PSOL. for(long iwfs=0; iwfs<aster->nwfs; iwfs++){ dadd(&mpsol->p[iwfs], 1, mreal, 1); } } pmerrm=0; if(istep>=parms->skyc.evlstart){/*performance evaluation*/ double res_ngs=dwdot(merr->p,parms->maos.mcc,merr->p); if(res_ngs>ngsol*100){ dfree(res); res=NULL; break; } { res->p[0]+=res_ngs; res->p[1]+=dwdot2(merr->p,parms->maos.mcc_tt,merr->p); double dot_oa=dwdot(merr->p, parms->maos.mcc_oa, merr->p); double dot_res_ideal=dwdot(merr->p, parms->maos.mcc_oa, PCOL(mideal,istep)); double dot_res_oa=0; for(int imod=0; imod<nmod; imod++){ dot_res_oa+=merr->p[imod]*IND(mideal_oa,imod,istep); } res->p[2]+=dot_oa-2*dot_res_ideal+2*dot_res_oa; res->p[4]+=dot_oa; } { double dot_oa_tt=dwdot2(merr->p, parms->maos.mcc_oa_tt, merr->p); /*Notice that mcc_oa_tt2 is 2x5 marix. */ double dot_res_ideal_tt=dwdot(merr->p, parms->maos.mcc_oa_tt2, PCOL(mideal,istep)); double dot_res_oa_tt=0; for(int imod=0; imod<2; imod++){ dot_res_oa_tt+=merr->p[imod]*IND(mideal_oa,imod,istep); } res->p[3]+=dot_oa_tt-2*dot_res_ideal_tt+2*dot_res_oa_tt; res->p[5]+=dot_oa_tt; } }//if evl if(istep<phystart || phystart<0){ /*Ztilt, noise free simulation for acquisition. */ dmm(&zgradc->m, 1, aster->gm, merr, "nn", 1);/*grad due to residual NGS mode. */ for(int iwfs=0; iwfs<aster->nwfs; iwfs++){ const int ipowfs=aster->wfs[iwfs].ipowfs; const long ng=parms->maos.nsa[ipowfs]*2; for(long ig=0; ig<ng; ig++){ zgradc->p[iwfs]->p[ig]+=aster->wfs[iwfs].ztiltout->p[istep*ng+ig]; } } for(int iwfs=0; iwfs<aster->nwfs; iwfs++){ int dtrati=(multirate?(int)dtrats->p[iwfs]:dtratc); if((istep+1) % dtrati==0){ dadd(&gradout->p[iwfs], 0, zgradc->p[iwfs], 1./dtrati); dzero(zgradc->p[iwfs]); if(noisy){ int idtrati=(multirate?(int)aster->idtrats->p[iwfs]:idtratc); dmat *nea=aster->wfs[iwfs].pistat->sanea->p[idtrati]; for(int i=0; i<nea->nx; i++){ gradout->p[iwfs]->p[i]+=nea->p[i]*randn(&aster->rand); } } pmerrm=merrm;//record output. } } }else{ /*Accumulate PSF intensities*/ for(long iwfs=0; iwfs<aster->nwfs; iwfs++){ const double thetax=aster->wfs[iwfs].thetax; const double thetay=aster->wfs[iwfs].thetay; const int ipowfs=aster->wfs[iwfs].ipowfs; const long nsa=parms->maos.nsa[ipowfs]; ccell* wvfout=aster->wfs[iwfs].wvfout[istep]; for(long iwvl=0; iwvl<nwvl; iwvl++){ double wvl=parms->maos.wvl[iwvl]; for(long isa=0; isa<nsa; isa++){ ccp(&wvfc->p[iwfs], IND(wvfout,isa,iwvl)); /*Apply NGS mode error to PSF. */ ngsmod2wvf(wvfc->p[iwfs], wvl, merr, powfs+ipowfs, isa, thetax, thetay, parms); cembedc(wvf->p[iwfs],wvfc->p[iwfs],0,C_FULL); cfft2(wvf->p[iwfs],-1); /*peak in corner. */ cabs22d(&psf[iwfs]->p[isa+nsa*iwvl], 1., wvf->p[iwfs], 1.); }/*isa */ }/*iwvl */ }/*iwfs */ /*Form detector image from accumulated PSFs*/ double igrad[2]; for(long iwfs=0; iwfs<aster->nwfs; iwfs++){ int dtrati=dtratc, idtrat=idtratc; if(multirate){//multirate idtrat=aster->idtrats->p[iwfs]; dtrati=dtrats->p[iwfs]; } if((istep+1) % dtrati == 0){/*has output */ dcellzero(ints[iwfs]); const int ipowfs=aster->wfs[iwfs].ipowfs; const long nsa=parms->maos.nsa[ipowfs]; for(long isa=0; isa<nsa; isa++){ for(long iwvl=0; iwvl<nwvl; iwvl++){ double siglev=aster->wfs[iwfs].siglev->p[iwvl]; ccpd(&otf->p[iwfs],psf[iwfs]->p[isa+nsa*iwvl]); cfft2i(otf->p[iwfs], 1); /*turn to OTF, peak in corner */ ccwm(otf->p[iwfs], powfs[ipowfs].dtf[iwvl].nominal); cfft2(otf->p[iwfs], -1); dspmulcreal(ints[iwfs]->p[isa]->p, powfs[ipowfs].dtf[iwvl].si, otf->p[iwfs]->p, siglev); } /*Add noise and apply matched filter. */ #if _OPENMP >= 200805 #pragma omp critical #endif switch(noisy){ case 0:/*no noise at all. */ break; case 1:/*both poisson and read out noise. */ { double bkgrnd=aster->wfs[iwfs].bkgrnd*dtrati; addnoise(ints[iwfs]->p[isa], &aster->rand, bkgrnd, bkgrnd, 0,0,IND(rnefs,idtrat,ipowfs)); } break; case 2:/*there is still poisson noise. */ addnoise(ints[iwfs]->p[isa], &aster->rand, 0, 0, 0,0,0); break; default: error("Invalid noisy\n"); } igrad[0]=0; igrad[1]=0; double pixtheta=parms->skyc.pixtheta[ipowfs]; if(parms->skyc.mtch){ dmulvec(igrad, mtche[iwfs]->p[isa], ints[iwfs]->p[isa]->p, 1); } if(!parms->skyc.mtch || fabs(igrad[0])>pixtheta || fabs(igrad[1])>pixtheta){ if(!parms->skyc.mtch){ warning2("fall back to cog\n"); }else{ warning_once("mtch is out of range\n"); } dcog(igrad, ints[iwfs]->p[isa], 0, 0, 0, 3*IND(rnefs,idtrat,ipowfs), 0); igrad[0]*=pixtheta; igrad[1]*=pixtheta; } gradout->p[iwfs]->p[isa]=igrad[0]; gradout->p[iwfs]->p[isa+nsa]=igrad[1]; }/*isa */ pmerrm=merrm; dcellzero(psf[iwfs]);/*reset accumulation.*/ }/*if iwfs has output*/ }/*for wfs*/ }/*if phystart */ //output to mreal after using it to ensure two cycle delay. if(st2t){//Type I or II control. if(st2t->mint->p[0]){//has output. dcp(&mreal, st2t->mint->p[0]->p[0]); } }else{//LQG control kalman_output(kalman, &mreal, 0, 1); } if(kalman){//LQG control int indk=0; //Form PSOL grads and obtain index to LQG M for(int iwfs=0; iwfs<aster->nwfs; iwfs++){ int dtrati=(multirate?(int)dtrats->p[iwfs]:dtratc); if((istep+1) % dtrati==0){ indk|=1<<iwfs; dmm(&gradout->p[iwfs], 1, aster->g->p[iwfs], mpsol->p[iwfs], "nn", 1./dtrati); dzero(mpsol->p[iwfs]); } } if(indk){ kalman_update(kalman, gradout->m, indk-1); } }else if(st2t){ if(pmerrm){ dmm(&merrm->p[0], 0, pgm, gradout->m, "nn", 1); } servo_filter(st2t, pmerrm);//do even if merrm is zero. to simulate additional latency } if(parms->skyc.dbg){ memcpy(PCOL(gradsave, istep), gradout->m->p, sizeof(double)*gradsave->nx); } }/*istep; */ } if(parms->skyc.dbg){ int dtrati=(multirate?(int)dtrats->p[0]:dtratc); writebin(gradsave,"%s/skysim_grads_aster%d_dtrat%d",dirsetup, aster->iaster,dtrati); writebin(mres,"%s/skysim_sim_mres_aster%d_dtrat%d",dirsetup,aster->iaster,dtrati); } dfree(mreal); dcellfree(mpsol); dfree(merr); dcellfree(merrm); dcellfree(zgradc); dcellfree(gradout); dfree(gradsave); if(hasphy){ dcellfreearr(psf, aster->nwfs); dcellfreearr(ints, aster->nwfs); ccellfree(wvf); ccellfree(wvfc); ccellfree(otf); free(mtche); } servo_free(st2t); /*dfree(mres); */ if(mresout) { *mresout=mres; }else{ dfree(mres); } dscale(res, 1./((nstep-parms->skyc.evlstart)*parms->skyc.navg)); return res; }