//--------------------------------------------------------- DMat& NDG2D::Lift2D() //--------------------------------------------------------- { // function [LIFT] = Lift2D() // Purpose : Compute surface to volume lift term for DG formulation DMat V1D,massEdge1,massEdge2,massEdge3; DVec faceR,faceS; Index1D J1(1,Nfp), J2(Nfp+1,2*Nfp), J3(2*Nfp+1,3*Nfp); DMat Emat(Np, Nfaces*Nfp); // face 1 faceR = r(Fmask(All,1)); V1D = Vandermonde1D(N, faceR); massEdge1 = inv(V1D*trans(V1D)); Emat(Fmask(All,1), J1) = massEdge1; // face 2 faceR = r(Fmask(All,2)); V1D = Vandermonde1D(N, faceR); massEdge2 = inv(V1D*trans(V1D)); Emat(Fmask(All,2), J2) = massEdge2; // face 3 faceS = s(Fmask(All,3)); V1D = Vandermonde1D(N, faceS); massEdge3 = inv(V1D*trans(V1D)); Emat(Fmask(All,3), J3) = massEdge3; // inv(mass matrix)*\I_n (L_i,L_j)_{edge_n} LIFT = V*(trans(V)*Emat); return LIFT; }
//--------------------------------------------------------- void CurvedINS2D::INSLiftDrag2D(double ra) //--------------------------------------------------------- { // function [Cd, Cl, dP, sw, stri] = INSLiftDrag2D(Ux, Uy, PR, ra, nu, time, tstep, Nsteps) // Purpose: compute coefficients of lift, drag and pressure drop at cylinder static FILE* fid; static DVec sw1, sw2; static int Nc=0, stri1=0, stri2=0; if (1 == tstep) { char buf[50]; sprintf(buf, "liftdraghistory%d.dat", N); fid = fopen(buf, "w"); fprintf(fid, "timeCdCldP = [...\n"); // Sample location and weights for pressure drop // Note: the VolkerCylinder test assumes the 2 // sample points (-ra, 0), (ra, 0), are vertices // located on the internal cylinder boundary Sample2D(-ra, 0.0, sw1, stri1); Sample2D( ra, 0.0, sw2, stri2); Nc = mapC.size()/Nfp; } bool bDo=false; if (1==tstep || tstep==Nsteps || !umMOD(tstep,10)) { bDo=true; } else if (time > 3.90 && time < 3.96) { bDo=true; } // catch C_d (3.9362) else if (time > 5.65 && time < 5.73) { bDo=true; } // catch C_l (5.6925) else if (time > 7.999) { bDo=true; } // catch dP (8.0) if (!bDo) return; DVec PRC, nxC, nyC, wv, tv; DMat dUxdx,dUxdy, dUydx,dUydy, MM1D, V1D; DMat dUxdxC,dUxdyC, dUydxC,dUydyC, hforce, vforce, sJC; double dP=0.0, Cd=0.0, Cl=0.0; dP = sw1.inner(PR(All,stri1)) - sw2.inner(PR(All,stri2)); // compute derivatives Grad2D(Ux, dUxdx,dUxdy); dUxdxC=dUxdx(vmapC); dUxdyC=dUxdy(vmapC); Grad2D(Uy, dUydx,dUydy); dUydxC=dUydx(vmapC); dUydyC=dUydy(vmapC); PRC=PR(vmapC); nxC=nx(mapC); nyC=ny(mapC); sJC=sJ(mapC); hforce = -PRC.dm(nxC) + nu*( nxC.dm( 2.0 *dUxdxC) + nyC.dm(dUydxC+dUxdyC) ); vforce = -PRC.dm(nyC) + nu*( nxC.dm(dUydxC+dUxdyC) + nyC.dm( 2.0 *dUydyC) ); hforce.reshape(Nfp, Nc); vforce.reshape(Nfp, Nc); sJC .reshape(Nfp, Nc); // compute weights for integrating (1,hforce) and (1,vforce) V1D = Vandermonde1D(N, r(Fmask(All,1))); MM1D = trans(inv(V1D)) / V1D; wv = MM1D.col_sums(); tv = wv*(sJC.dm(hforce)); Cd=tv.sum(); // Compute drag coefficient tv = wv*(sJC.dm(vforce)); Cl=tv.sum(); // Compute lift coefficient // Output answers for plotting fprintf(fid, "%15.14f %15.14f %15.14f %15.14f;...\n", time, Cd/ra, Cl/ra, dP); fflush(fid); // LOG report if (1==tstep || tstep==Nsteps || !umMOD(tstep,Nreport)) { umLOG(1, "%7d %6.3lf %9.5lf %10.6lf %9.5lf\n", tstep, time, Cd/ra, Cl/ra, dP); } if (tstep==Nsteps) { fprintf(fid, "];\n"); fclose(fid); fid=NULL; } }
//--------------------------------------------------------- void NDG2D::MakeCylinder2D ( const IMat& faces, double ra, double xo, double yo ) //--------------------------------------------------------- { // Function: MakeCylinder2D(faces, ra, xo, yo) // Purpose: Use Gordon-Hall blending with an isoparametric // map to modify a list of faces so they conform // to a cylinder of radius r centered at (xo,yo) int NCurveFaces = faces.num_rows(); IVec vflag(VX.size()); IMat VFlag(EToV.num_rows(),EToV.num_cols()); int n=0,k=0,f=0,v1=0,v2=0; double theta1=0.0,theta2=0.0,x1=0.0,x2=0.0,y1=0.0,y2=0.0; double newx1=0.0,newx2=0.0,newy1=0.0,newy2=0.0; DVec vr, fr, theta, fdx,fdy, vdx, vdy, blend,numer,denom; IVec va,vb,vc, ks, Fm_f, ids; DMat Vface, Vvol; for (n=1; n<=NCurveFaces; ++n) { // move vertices of faces to be curved onto circle k = faces(n,1); f = faces(n,2); v1 = EToV(k, f); v2 = EToV(k, umMOD(f,Nfaces)+1); theta1 = atan2(VY(v1),VX(v1)); theta2 = atan2(VY(v2),VX(v2)); newx1 = xo + ra*cos(theta1); newy1 = yo + ra*sin(theta1); newx2 = xo + ra*cos(theta2); newy2 = yo + ra*sin(theta2); // update mesh vertex locations VX(v1) = newx1; VX(v2) = newx2; VY(v1) = newy1; VY(v2) = newy2; // store modified vertex numbers vflag(v1) = 1; vflag(v2) = 1; } // map modified vertex flag to each element //vflag = vflag(EToV); VFlag.resize(EToV); VFlag.set_map(EToV, vflag); // (map, values) // locate elements with at least one modified vertex //ks = find(sum(vflag,2)>0); ks = find(VFlag.row_sums(), '>', 0); // build coordinates of all the corrected nodes IMat VA=EToV(ks,1), VB=EToV(ks,2), VC=EToV(ks,3); // FIXME: loading 2D mapped data into 1D vectors int Nr=VA.num_rows(); va.copy(Nr,VA.data()); vb.copy(Nr,VB.data()); vc.copy(Nr,VC.data()); // Note: outer products of (Vector,MappedRegion1D) x(All,ks) = 0.5*(-(r+s)*VX(va)+(1.0+r)*VX(vb)+(1.0+s)*VX(vc)); y(All,ks) = 0.5*(-(r+s)*VY(va)+(1.0+r)*VY(vb)+(1.0+s)*VY(vc)); // deform specified faces for (n=1; n<=NCurveFaces; ++n) { k = faces(n,1); f = faces(n,2); // find vertex locations for this face and tangential coordinate if (f==1) { v1=EToV(k,1); v2=EToV(k,2); vr=r; } else if (f==2) { v1=EToV(k,2); v2=EToV(k,3); vr=s; } else if (f==3) { v1=EToV(k,1); v2=EToV(k,3); vr=s; } fr = vr(Fmask(All,f)); x1 = VX(v1); y1 = VY(v1); x2 = VX(v2); y2 = VY(v2); // move vertices at end points of this face to the cylinder theta1 = atan2(y1-yo, x1-xo); theta2 = atan2(y2-yo, x2-xo); // check to make sure they are in the same quadrant if ((theta2 > 0.0) && (theta1 < 0.0)) { theta1 += 2*pi; } if ((theta1 > 0.0) && (theta2 < 0.0)) { theta2 += 2*pi; } // distribute N+1 nodes by arc-length along edge theta = 0.5*theta1*(1.0-fr) + 0.5*theta2*(1.0+fr); // evaluate deformation of coordinates fdx = xo + ra*apply(cos,theta) - x(Fmask(All,f),k); fdy = yo + ra*apply(sin,theta) - y(Fmask(All,f),k); // build 1D Vandermonde matrix for face nodes and volume nodes Vface = Vandermonde1D(N, fr); Vvol = Vandermonde1D(N, vr); // compute unblended volume deformations vdx = Vvol * (Vface|fdx); vdy = Vvol * (Vface|fdy); // blend deformation and increment node coordinates ids = find(abs(1.0-vr), '>', 1e-7); // warp and blend denom = 1.0-vr(ids); if (1==f) { numer = -(r(ids)+s(ids)); } else if (2==f) { numer = (r(ids)+1.0 ); } else if (3==f) { numer = -(r(ids)+s(ids)); } blend = numer.dd(denom); x(ids,k) += (blend.dm(vdx(ids))); y(ids,k) += (blend.dm(vdy(ids))); } // repair other coordinate dependent information Fx = x(Fmask, All); Fy = y(Fmask, All); ::GeometricFactors2D(x,y,Dr,Ds, rx,sx,ry,sy,J); Normals2D(); Fscale = sJ.dd(J(Fmask,All)); }
//--------------------------------------------------------- DVec& NDG2D::PoissonIPDGbc2D (DVec& ubc, //[in] DVec& qbc //[in] ) //--------------------------------------------------------- { // function [OP] = PoissonIPDGbc2D() // Purpose: Set up the discrete Poisson matrix directly // using LDG. The operator is set up in the weak form // build DG derivative matrices int max_OP = (K*Np*Np*(1+Nfaces)); // initialize parameters DVec faceR("faceR"), faceS("faceS"); DMat V1D("V1D"), Dx("Dx"),Dy("Dy"), Dn1("Dn1"), mmE_Fm1("mmE(:,Fm1)"); IVec Fm("Fm"), Fm1("Fm1"), fidM("fidM"); double lnx=0.0,lny=0.0,lsJ=0.0,hinv=0.0,gtau=0.0; int i=0,k1=0,f1=0,id=0; IVec i1_Nfp = Range(1,Nfp); double N1N1 = double((N+1)*(N+1)); // build local face matrices DMat massEdge[4]; // = zeros(Np,Np,Nfaces); for (i=1; i<=Nfaces; ++i) { massEdge[i].resize(Np,Np); } // face mass matrix 1 Fm = Fmask(All,1); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[1](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 2 Fm = Fmask(All,2); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[2](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 3 Fm = Fmask(All,3); faceS = s(Fm); V1D = Vandermonde1D(N, faceS); massEdge[3](Fm,Fm) = inv(V1D*trans(V1D)); // build DG right hand side DVec* pBC = new DVec(Np*K, "bc", OBJ_temp); DVec& bc = (*pBC); // reference, for syntax //////////////////////////////////////////////////////////////// umMSG(1, "\n ==> {OP} assembly [bc]: "); for (k1=1; k1<=K; ++k1) { if (! (k1%100)) { umMSG(1, "%d, ",k1); } // rows1 = outer(Range((k1-1)*Np+1,k1*Np), Ones(NGauss)); // Build element-to-element parts of operator for (f1=1; f1<=Nfaces; ++f1) { if (BCType(k1,f1)) { ////////////////////////added by Kevin /////////////////////////////// Fm1 = Fmask(All,f1); fidM = (k1-1)*Nfp*Nfaces + (f1-1)*Nfp + i1_Nfp; id = 1+(f1-1)*Nfp + (k1-1)*Nfp*Nfaces; lnx = nx(id); lny = ny(id); lsJ = sJ(id); hinv = Fscale(id); Dx = rx(1,k1)*Dr + sx(1,k1)*Ds; Dy = ry(1,k1)*Dr + sy(1,k1)*Ds; Dn1 = lnx*Dx + lny*Dy; //mmE = lsJ*massEdge(:,:,f1); //bc(All,k1) += (gtau*mmE(All,Fm1) - Dn1'*mmE(All,Fm1))*ubc(fidM); mmE_Fm1 = massEdge[f1](All,Fm1); mmE_Fm1 *= lsJ; gtau = 10*N1N1*hinv; // set penalty scaling //bc(All,k1) += (gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1) * ubc(fidM); switch(BCType(k1,f1)){ case BC_Dirichlet: bc(Np*(k1-1)+Range(1,Np)) += (gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1)*ubc(fidM); break; case BC_Neuman: bc(Np*(k1-1)+Range(1,Np)) += mmE_Fm1*qbc(fidM); break; default: std::cout<<"warning: boundary condition is incorrect"<<std::endl; } } } } return bc; }
void NDG2D::PoissonIPDGbc2D( CSd& spOP //[out] sparse operator ) { // function [OP] = PoissonIPDGbc2D() // Purpose: Set up the discrete Poisson matrix directly // using LDG. The operator is set up in the weak form // build DG derivative matrices int max_OP = (K*Np*Np*(1+Nfaces)); //initialize parameters DVec faceR("faceR"), faceS("faceS"); IVec Fm("Fm"), Fm1("Fm1"), fidM("fidM"); DMat V1D("V1D"); int i=0; // build local face matrices DMat massEdge[4]; // = zeros(Np,Np,Nfaces); for (i=1; i<=Nfaces; ++i) { massEdge[i].resize(Np,Np); } // face mass matrix 1 Fm = Fmask(All,1); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[1](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 2 Fm = Fmask(All,2); faceR = r(Fm); V1D = Vandermonde1D(N, faceR); massEdge[2](Fm,Fm) = inv(V1D*trans(V1D)); // face mass matrix 3 Fm = Fmask(All,3); faceS = s(Fm); V1D = Vandermonde1D(N, faceS); massEdge[3](Fm,Fm) = inv(V1D*trans(V1D)); //continue initialize parameters DMat Dx("Dx"),Dy("Dy"), Dn1("Dn1"), mmE_Fm1("mmE(:,Fm1)"); double lnx=0.0,lny=0.0,lsJ=0.0,hinv=0.0,gtau=0.0; int k1=0,f1=0,id=0; IVec i1_Nfp = Range(1,Nfp); double N1N1 = double((N+1)*(N+1)); // "OP" triplets (i,j,x), extracted to {Ai,Aj,Ax} IVec OPi(max_OP),OPj(max_OP), Ai,Aj; DVec OPx(max_OP), Ax; IMat rows1, cols1; Index1D entries; DMat OP11(Np,Nfp, 0.0); // global node numbering entries.reset(1,Np*Nfp); cols1 = outer(Ones(Np), Range(1,Nfp)); umMSG(1, "\n ==> {OP} assembly [bc]: "); for (k1=1; k1<=K; ++k1) { if (! (k1%100)) { umMSG(1, "%d, ",k1); } rows1 = outer(Range((k1-1)*Np+1,k1*Np), Ones(Nfp)); // Build element-to-element parts of operator for (f1=1; f1<=Nfaces; ++f1) { if (BCType(k1,f1)) { ////////////////////////added by Kevin /////////////////////////////// Fm1 = Fmask(All,f1); fidM = (k1-1)*Nfp*Nfaces + (f1-1)*Nfp + i1_Nfp; id = 1+(f1-1)*Nfp + (k1-1)*Nfp*Nfaces; lnx = nx(id); lny = ny(id); lsJ = sJ(id); hinv = Fscale(id); Dx = rx(1,k1)*Dr + sx(1,k1)*Ds; Dy = ry(1,k1)*Dr + sy(1,k1)*Ds; Dn1 = lnx*Dx + lny*Dy; //mmE = lsJ*massEdge(:,:,f1); //bc(All,k1) += (gtau*mmE(All,Fm1) - Dn1'*mmE(All,Fm1))*ubc(fidM); mmE_Fm1 = massEdge[f1](All,Fm1); mmE_Fm1 *= lsJ; gtau = 10*N1N1*hinv; // set penalty scaling //bc(All,k1) += (gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1) * ubc(fidM); switch(BCType(k1,f1)){ case BC_Dirichlet: OP11 = gtau*mmE_Fm1 - trans(Dn1)*mmE_Fm1; break; case BC_Neuman: OP11 = mmE_Fm1; break; default: std::cout<<"warning: boundary condition is incorrect"<<std::endl; } OPi(entries)=rows1; OPj(entries)=cols1; OPx(entries)=OP11; entries += (Np*Nfp); } cols1 += Nfp; } } umMSG(1, "\n ==> {OPbc} to sparse\n"); entries.reset(1, entries.hi()-(Np*Nfp)); // extract triplets from large buffers Ai=OPi(entries); Aj=OPj(entries); Ax=OPx(entries); // These arrays can be HUGE, so force deallocation OPi.Free(); OPj.Free(); OPx.Free(); // return 0-based sparse result Ai -= 1; Aj -= 1; //------------------------------------------------------- // This operator is not symmetric, and will NOT be // factorised, only used to create reference RHS's: // // refrhsbcPR = spOP1 * bcPR; // refrhsbcUx = spOP2 * bcUx; // refrhsbcUy = spOP2 * bcUy; // // Load ALL elements (both upper and lower triangles): //------------------------------------------------------- spOP.load(Np*K, Nfp*Nfaces*K, Ai,Aj,Ax, sp_All,false, 1e-15,true); Ai.Free(); Aj.Free(); Ax.Free(); umMSG(1, " ==> {OPbc} ready.\n"); #if (1) // check on original estimates for nnx umMSG(1, " ==> max_OP: %12d\n", max_OP); umMSG(1, " ==> nnz_OP: %12d\n", entries.hi()); #endif }
startUp1d::startUp1d(int N,int K) { NODETOL=1e-10; r = JacobiGL(0,0,N); Np = N+1; Nfp = 1; Nfaces =2; V = Vandermonde1D(N,r); Dr = Dmatrix1D(N,r,V); invV = V.inverse(); MatrixXd Emat(Np,Nfaces*Nfp); Emat.setZero(Np,2); Emat(0,0)=1.0; Emat(Np-1,1)=1.0; LIFT = V*V.transpose()*Emat; EToV.setZero(K,2); VX.setZero(K+1); EToV = meshGen(0.0,1.0,K,Nv,VX); ArrayXi va,vb; va = EToV.col(0); vb=EToV.col(1); VectorXd v1(K),v2(K); { int temp=0; for(int i=0; i<K; i++) { temp = va(i); v1(i)= VX(temp-1); temp = vb(i); v2(i)= VX(temp-1); } } VectorXd temp1; temp1.setOnes(N+1,1); x =temp1*v1.transpose()+ 0.5*((r.array()+1).matrix())*(v2-v1).transpose(); J = Dr*x; rx = 1.0/J.array(); Fmask.setZero(2); Fmask(0)=0; Fmask(1)=Np-1; std::cout<<"Fx"<<Fmask<<"\n"; Fx.setZero(2,K); Fx.row(0)=x.row(0); Fx.row(1)=x.row(Np-1); std::cout<<"Fx"<<Fx<<"\n"; nx=Normals1D(Nfp,Nfaces,K); Fscale.setZero(2,K); Fscale.row(0)=J.row(0); Fscale.row(1)=J.row(Np-1); Fscale = 1.0/Fscale.array(); EToEF = Connect1DEToE(EToV); std::cout<<"Fscale "<<"\n"<<Fscale<<"\n"; std::cout<<"Fscale "<<"\n"<<EToEF<<"\n"; }