// some unit tests of the macromesh code int TestMacroMesh(void) { MacroMesh m; int param[]={4, 4, 4, 1, 1, 1, 0}; // test gmsh file reading ReadMacroMesh(&m, "test/testmacromesh.msh"); BuildConnectivity(&m); CheckMacroMesh(&m, param); PrintMacroMesh(&m); int test = (m.nbelems == 5); test = (test && m.nbnodes == 50); // test search methods real xphy[3]={1,1.1,0.5}; real xref[3]; test= test && IsInElem(&m,0,xphy,xref); printf("xphy=%f %f %f xref=%f %f %f \n",xphy[0],xphy[1],xphy[2], xref[0],xref[1],xref[2]); xphy[2]=-0.5; test= test && !IsInElem(&m,0,xphy,xref); int num=NumElemFromPoint(&m,xphy,NULL); printf("xphy=%f %f %f is in elem=%d\n",xphy[0],xphy[1],xphy[2],num); test=test && (num == -1); xphy[2]=0.5; num=NumElemFromPoint(&m,xphy,NULL); printf("xphy=%f %f %f is in elem=%d\n",xphy[0],xphy[1],xphy[2],num); test=test && (num == 0); real xphy2[3]={1,0,0.33}; num=NumElemFromPoint(&m,xphy2,NULL); printf("xphy=%f %f %f is in elem=%d\n",xphy2[0],xphy2[1],xphy2[2],num); test=test && (num == 3); return test; }
void CreateCoil2DParticles(PIC* pic,MacroMesh *m){ real delta=1; // coil radius real current=1; // electric current real v0=1; // tangential particle velocity real pi=4*atan(1.); pic->weight= pi*delta*current/v0/pic->nbparts; for(int np=0;np<pic->nbparts;np++){ real xphi[3]; real rax=corput(np+1,5,3); xphi[0]=delta*cos(rax*2*pi); xphi[1]=delta*sin(rax*2*pi); xphi[2]=0.5; real vx=-v0*sin(rax*2*pi); real vy=v0*cos(rax*2*pi); real vz=0; real weight=2*pi*delta*current/v0/pic->nbparts; real xref[3]; int num_elem=NumElemFromPoint(m,xphi,xref); pic->cell_id[np]=num_elem; pic->old_cell_id[np]=num_elem; pic->xv[np*6+0]=xref[0]; pic->xv[np*6+1]=xref[1]; pic->xv[np*6+2]=xref[2]; pic->xv[np*6+3]=vx; pic->xv[np*6+4]=vy; pic->xv[np*6+5]=vz; } }
void PushParticles(field *f,PIC* pic){ for(int i=0;i<pic->nbparts;i++) { // jacobian of tau at the particle real physnode[20][3]; int ie=pic->cell_id[i]; if (ie >=0) { real w[f->model.m]; InterpField(f,pic->cell_id[i],&(pic->xv[6*i]),w); for(int inoloc = 0; inoloc < 20; inoloc++) { int ino = f->macromesh.elem2node[20*ie+inoloc]; physnode[inoloc][0] = f->macromesh.node[3 * ino + 0]; physnode[inoloc][1] = f->macromesh.node[3 * ino + 1]; physnode[inoloc][2] = f->macromesh.node[3 * ino + 2]; } real dtau[3][3],codtau[3][3]; real xphy[3]; Ref2Phy(physnode, // phys. nodes &(pic->xv[6*i]), // xref NULL, -1, // dpsiref, ifa xphy, dtau, // xphy, dtau codtau, NULL, NULL); // codtau, dpsi, vnds real det = dot_product(dtau[0], codtau[0]); //printf("w=%f %f %f %f\n",w[0],w[1],w[2],w[3]); // 2D motion real vref[3]; pic->xv[6*i+3+0] +=pic->dt * (w[0]+w[2]*pic->xv[6*i+4]); pic->xv[6*i+3+1] +=pic->dt * (w[1]-w[2]*pic->xv[6*i+3]); pic->xv[6*i+3+2] +=0; for(int ii=0;ii<3;ii++){ vref[ii]=0; for(int jj=0;jj<3;jj++){ vref[ii] += codtau[jj][ii] * pic->xv[6*i+3+jj] / det ; } } // TO DO: check if the particle is in a new element pic->xv[6*i+0]+=pic->dt * vref[0]; pic->xv[6*i+1]+=pic->dt * vref[1]; pic->xv[6*i+2]+=pic->dt * vref[2]; bool is_out = (pic->xv[6*i+0] < 0 || pic->xv[6*i+0] > 1) || (pic->xv[6*i+1] < 0 || pic->xv[6*i+1] > 1) || (pic->xv[6*i+2] < 0 || pic->xv[6*i+2] > 1); //is_out = false; if (is_out) { Ref2Phy(physnode, // phys. nodes &(pic->xv[6*i]), // xref NULL, -1, // dpsiref, ifa xphy, NULL, // xphy, dtau NULL, NULL, NULL); // codtau, dpsi, vnds int old=pic->cell_id[i]; printf("oldref=%f %f %f \n",pic->xv[6*i+0], pic->xv[6*i+1],pic->xv[6*i+2]); pic->cell_id[i]= NumElemFromPoint(&(f->macromesh), xphy, &(pic->xv[6*i])); if (pic->cell_id[i] != -1) pic->old_cell_id[i]=pic->cell_id[i]; printf("newref=%f %f %f \n",pic->xv[6*i+0], pic->xv[6*i+1],pic->xv[6*i+2]); printf("change elem: %d -> %d \n",old,pic->cell_id[i]); } } // if ie >= 0 } }
void CreateParticles(PIC* pic,MacroMesh *m){ int k1[7]={3,5,7,11,13,17,19}; int k2[7]={2,3,5,7,11,13,17}; int n=0; int np=0; real vt=1; real xp[3]; real vp[3]; //srand(time(NULL)); //int r = rand(); while(np < pic->nbparts){ for(int idim=0;idim<3;idim++){ real r=corput(n,k1[idim],k2[idim]); //real r=rand(); //r/=RAND_MAX; xp[idim]=m->xmin[idim]+r* (m->xmax[idim]-m->xmin[idim]); } real xref[3]; int num_elem=NumElemFromPoint(m,xp,xref); //printf("numelem=%d %f %f %f\n",num_elem,xp[0],xp[1],xp[2]); if (num_elem != -1) { /* pic->xv[np*6+0]=xp[0]; */ /* pic->xv[np*6+1]=xp[1]; */ /* pic->xv[np*6+2]=xp[2]; */ pic->xv[np*6+0]=xref[0]; pic->xv[np*6+1]=xref[1]; pic->xv[np*6+2]=xref[2]; real vp[3]={1,0,0}; BoxMuller3d(vp,k1+3,k2+3); /* for(int idim=0;idim<3;idim++){ vp[idim]=vt*(corput(n,k1[idim+3],k2[idim+3])-0.5); vp[idim]=0; //!!!!!!!!!!!!!!!!!! }*/ pic->xv[np*6+3]=vp[0]; pic->xv[np*6+4]=vp[1]; pic->xv[np*6+5]=vp[2]; pic->cell_id[np]=num_elem; pic->old_cell_id[np]=num_elem; np++; } n++; } }
// some unit tests of the macromesh code int TestPICAccumulate(void) { MacroMesh m; bool test=true; int param[]={4, 4, 4, 1, 1, 1, 0}; field f; init_empty_field(&f); // test gmsh file reading ReadMacroMesh(&(f.macromesh), "test/testmacromesh.msh"); BuildConnectivity(&(f.macromesh)); CheckMacroMesh(&(f.macromesh), param); //PrintMacroMesh(&m); PIC pic; InitPIC(&pic,1); CreateParticles(&pic,&(f.macromesh)); PlotParticles(&pic,&(f.macromesh)); f.model.m = 7; // num of conservative variables /* f.model.NumFlux = Maxwell2DNumFlux; */ /* f.model.BoundaryFlux = Maxwell2DBoundaryFlux; */ f.model.InitData = Maxwell2DConstInitData; /* f.model.ImposedData = Maxwell2DImposedData; */ f.varindex = GenericVarindex; f.interp.interp_param[0] = f.model.m; f.interp.interp_param[1] = 1; // x direction degree f.interp.interp_param[2] = 1; // y direction degree f.interp.interp_param[3] = 0; // z direction degree f.interp.interp_param[4] = 1; // x direction refinement f.interp.interp_param[5] = 1; // y direction refinement f.interp.interp_param[6] = 1; // z direction refinement Initfield(&f); // place the particle at (0,1,0) and v=(1,0,0) pic.xv[0]=0; pic.xv[1]=0; pic.xv[2]=0.5; real xref[3]; pic.cell_id[0]=NumElemFromPoint(&f.macromesh,pic.xv,xref); pic.xv[0]=xref[0]; pic.xv[1]=xref[1]; pic.xv[2]=xref[2]; pic.xv[3]=1; pic.xv[4]=0; pic.xv[5]=0; PlotParticles(&pic,&(f.macromesh)); int ie=2; int ipg=2; int iv=4; int imem=f.varindex(f.interp_param, ie, ipg, iv); AccumulateParticles(&pic,&f); printf("w=%f wex=%f\n",f.wn[imem],1/1.96); test = test && (fabs(f.wn[imem]-1/1.96) < 1e-8); Displayfield(&f); return test; }
// Build other connectivity arrays void BuildConnectivity(MacroMesh* m) { printf("Build connectivity...\n"); real *bounds = malloc(6 * sizeof(real)); macromesh_bounds(m, bounds); printf("bounds: %f, %f, %f, %f, %f, %f\n", bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]); printf("bounds: %f, %f, %f, %f, %f, %f\n", m->xmin[0],m->xmax[0], m->xmin[1],m->xmax[1], m->xmin[2],m->xmax[2] ); // Build a list of faces each face is made of four corners of the // hexaedron mesh Face4Sort *face = malloc(6 * sizeof(Face4Sort) * m->nbelems); build_face(m, face); build_elem2elem(m, face); free(face); build_node2elem(m); // check /* for(int ie = 0;ie<m->nbelems;ie++) { */ /* for(int ifa = 0;ifa<6;ifa++) { */ /* printf("elem=%d face=%d, voisin=%d\n", */ /* ie,ifa,m->elem2elem[ifa+6*ie]); */ /* } */ /* } */ if(m->is2d) suppress_zfaces(m); if(m->is1d) { suppress_zfaces(m); suppress_yfaces(m); } // update connectivity if the mesh is periodic // in some directions real diag[3][3]={1,0,0, 0,1,0, 0,0,1}; for (int ie = 0; ie < m->nbelems; ie++) { real physnode[20][3]; for(int inoloc = 0; inoloc < 20; inoloc++) { int ino = m->elem2node[20 * ie + inoloc]; physnode[inoloc][0] = m->node[3 * ino + 0]; physnode[inoloc][1] = m->node[3 * ino + 1]; physnode[inoloc][2] = m->node[3 * ino + 2]; } for(int ifa = 0; ifa < 6; ifa++) { if (m->elem2elem[6 * ie + ifa] < 0){ real xpgref[3],xpgref_in[3]; int ipgf=0; int param2[7]={0,0,0,1,1,1,0}; ref_pg_face(param2, ifa, ipgf, xpgref, NULL, xpgref_in); real dtau[3][3],xpg_in[3]; real codtau[3][3],vnds[3]={0,0,0}; Ref2Phy(physnode, xpgref_in, NULL, ifa, // dpsiref,ifa xpg_in, dtau, codtau, NULL, vnds); // codtau,dpsi,vnds Normalize(vnds); vnds[0]=fabs(vnds[0]); vnds[1]=fabs(vnds[1]); vnds[2]=fabs(vnds[2]); int dim=0; while(Dist(vnds,diag[dim]) > 1e-2 && dim<3) dim++; //assert(dim < 3); //printf("xpg_in_before=%f\n",xpg_in[dim]); if (dim < 3 && m->period[dim] > 0){ //if (xpg_in[dim] > m->period[dim]){ if (xpg_in[dim] > m->xmax[dim]){ xpg_in[dim] -= m->period[dim]; //printf("xpg_in_after=%f\n",xpg_in[dim]); } //else if (xpg_in[dim] < 0){ else if (xpg_in[dim] < m->xmin[dim]){ xpg_in[dim] += m->period[dim]; //printf("xpg_in_after=%f\n",xpg_in[dim]); } else { //printf("xpg_in=%f\n",xpg_in[dim]); assert(1==2); } m->elem2elem[6 * ie + ifa] = NumElemFromPoint(m,xpg_in,NULL); /* printf("ie=%d ifa=%d numelem=%d vnds=%f %f %f xpg_in=%f %f %f \n", */ /* ie,ifa,NumElemFromPoint(m,xpg_in,NULL), */ /* vnds[0],vnds[1],vnds[2], */ /* xpg_in[0],xpg_in[1],xpg_in[2]); */ } } } } // now, update the face2elem connectivity (because elem2elem has changed) for(int ifa = 0; ifa < m->nbfaces; ifa++) { int ieL = m->face2elem[4 * ifa + 0]; int locfaL = m->face2elem[4 * ifa + 1]; int ieR = m->face2elem[4 * ifa + 2]; int locfaR = m->face2elem[4 * ifa + 3]; int ieR2=m->elem2elem[6 * ieL + locfaL]; if (ieR != ieR2){ assert(ieR == -1); int opp[6]={2,3,0,1,5,4}; if (locfaL == 0 || locfaL == 1 || locfaL == 4) { m->face2elem[4 * ifa + 2] = ieR2; m->face2elem[4 * ifa + 3] = opp[locfaL]; } else { // mark the face for suppression m->face2elem[4 * ifa + 0] = -1; } } } suppress_double_faces(m); //assert(1==5); free(bounds); m->connec_ok = true; /* #ifdef _PERIOD */ /* assert(m->is1d); // TODO : generalize to 2D */ /* assert(m->nbelems==1); */ /* // faces 1 and 3 point to the same unique macrocell */ /* m->elem2elem[1+6*0]=0; */ /* m->elem2elem[3+6*0]=0; */ /* #endif */ }