// GAUSS-LOBATTO QUADRATURE ALONG EDGE void SetEdgeDataGL_Unst(const mesh& Mesh, int NumQuadPoints, int NumBasisOrder, edge_data_Unst* EdgeData) { // Quick error check if (NumQuadPoints<2 || NumQuadPoints>6 || NumBasisOrder<1 || NumBasisOrder>5) { printf(" \n"); printf(" Error in SetEdgeData_Unst.cpp \n"); printf(" NumQuadPoints must be 2,3,4,5, or 6.\n"); printf(" NumBasisOrder must be 1,2,3,4, or 5.\n"); printf(" NumQuadPoints = %i\n",NumQuadPoints); printf(" NumBasisOrder = %i\n",NumBasisOrder); printf("\n"); exit(1); } // --------------------------------- // Set quadrature weights and points // --------------------------------- switch( NumQuadPoints ) { case 2: EdgeData->GL_wgts1d->set(1, 1.0 ); EdgeData->GL_wgts1d->set(2, 1.0 ); EdgeData->GL_xpts1d->set(1, 1.0 ); EdgeData->GL_xpts1d->set(2, -1.0 ); break; case 3: EdgeData->GL_wgts1d->set(1, onethird ); EdgeData->GL_wgts1d->set(2, 4.0*onethird ); EdgeData->GL_wgts1d->set(3, onethird ); EdgeData->GL_xpts1d->set(1, 1.0 ); EdgeData->GL_xpts1d->set(2, 0.0 ); EdgeData->GL_xpts1d->set(3, -1.0 ); break; case 4: EdgeData->GL_wgts1d->set(1, 0.5*onethird ); EdgeData->GL_wgts1d->set(2, 2.5*onethird ); EdgeData->GL_wgts1d->set(3, 2.5*onethird ); EdgeData->GL_wgts1d->set(4, 0.5*onethird ); EdgeData->GL_xpts1d->set(1, 1.0 ); EdgeData->GL_xpts1d->set(2, osq5 ); EdgeData->GL_xpts1d->set(3, -osq5 ); EdgeData->GL_xpts1d->set(4, -1.0 ); break; case 5: EdgeData->GL_wgts1d->set(1, 0.1 ); EdgeData->GL_wgts1d->set(2, 49.0/90.0 ); EdgeData->GL_wgts1d->set(3, 32.0/45.0 ); EdgeData->GL_wgts1d->set(4, 49.0/90.0 ); EdgeData->GL_wgts1d->set(5, 0.1 ); EdgeData->GL_xpts1d->set(1, 1.0 ); EdgeData->GL_xpts1d->set(2, sq3*osq7 ); EdgeData->GL_xpts1d->set(3, 0.0 ); EdgeData->GL_xpts1d->set(4, -sq3*osq7 ); EdgeData->GL_xpts1d->set(5, -1.0 ); break; case 6: EdgeData->GL_wgts1d->set(1, 0.2*onethird ); EdgeData->GL_wgts1d->set(2, (1.4 - 0.1*sq7)*onethird ); EdgeData->GL_wgts1d->set(3, (1.4 + 0.1*sq7)*onethird ); EdgeData->GL_wgts1d->set(4, (1.4 + 0.1*sq7)*onethird ); EdgeData->GL_wgts1d->set(5, (1.4 - 0.1*sq7)*onethird ); EdgeData->GL_wgts1d->set(6, 0.2*onethird ); EdgeData->GL_xpts1d->set(1, 1.0 ); EdgeData->GL_xpts1d->set(2, (1/21.0)*sqrt(147.0+42.0*sq7) ); EdgeData->GL_xpts1d->set(3, (1/21.0)*sqrt(147.0-42.0*sq7) ); EdgeData->GL_xpts1d->set(4, -(1/21.0)*sqrt(147.0-42.0*sq7) ); EdgeData->GL_xpts1d->set(5, -(1/21.0)*sqrt(147.0+42.0*sq7) ); EdgeData->GL_xpts1d->set(6, -1.0 ); break; } // --------------------------------- // Legendre basis functions on the // left and right of each edge // --------------------------------- const int NumEdges = Mesh.get_NumEdges(); const int NumBasisComps = (NumBasisOrder*(NumBasisOrder+1))/2; dTensor1 xp1(3); dTensor1 yp1(3); dTensor1 xp2(3); dTensor1 yp2(3); dTensor1 xy1(2); dTensor1 xy2(2); dTensor1 mu1(NumBasisComps); dTensor1 mu2(NumBasisComps); for (int i=1; i<=NumEdges; i++) { // Get edge information const double x1 = Mesh.get_edge(i,1); const double y1 = Mesh.get_edge(i,2); const double x2 = Mesh.get_edge(i,3); const double y2 = Mesh.get_edge(i,4); const int e1 = Mesh.get_eelem(i,1); const int e2 = Mesh.get_eelem(i,2); // Get element information about // the two elements that meet at // the current edge const double Area1 = Mesh.get_area_prim(e1); const double Area2 = Mesh.get_area_prim(e2); for (int k=1; k<=3; k++) { xp1.set(k, Mesh.get_node(Mesh.get_tnode(e1,k),1) ); yp1.set(k, Mesh.get_node(Mesh.get_tnode(e1,k),2) ); xp2.set(k, Mesh.get_node(Mesh.get_tnode(e2,k),1) ); yp2.set(k, Mesh.get_node(Mesh.get_tnode(e2,k),2) ); } const double xc1 = (xp1.get(1) + xp1.get(2) + xp1.get(3))/3.0; const double yc1 = (yp1.get(1) + yp1.get(2) + yp1.get(3))/3.0; const double xc2 = (xp2.get(1) + xp2.get(2) + xp2.get(3))/3.0; const double yc2 = (yp2.get(1) + yp2.get(2) + yp2.get(3))/3.0; // quadrature points on the edge for (int m=1; m<=NumQuadPoints; m++) { // Take integration point s (in [-1,1]) // and map to physical domain const double s = EdgeData->GL_xpts1d->get(m); const double x = x1 + 0.5*(s+1.0)*(x2-x1); const double y = y1 + 0.5*(s+1.0)*(y2-y1); // Take physical point (x,y) // and map into the coordinates // of the two triangles that are // adjacent to the current edge xy1.set(1, ((yp1.get(3)-yp1.get(1))*(x-xc1) + (xp1.get(1)-xp1.get(3))*(y-yc1))/(2.0*Area1) ); xy1.set(2, ((yp1.get(1)-yp1.get(2))*(x-xc1) + (xp1.get(2)-xp1.get(1))*(y-yc1))/(2.0*Area1) ); xy2.set(1, ((yp2.get(3)-yp2.get(1))*(x-xc2) + (xp2.get(1)-xp2.get(3))*(y-yc2))/(2.0*Area2) ); xy2.set(2, ((yp2.get(1)-yp2.get(2))*(x-xc2) + (xp2.get(2)-xp2.get(1))*(y-yc2))/(2.0*Area2) ); // Evaluate monomials at locations xy1 double xi = xy1.get(1); double xi2 = xi*xi; double xi3 = xi*xi2; double xi4 = xi*xi3; double eta = xy1.get(2); double eta2 = eta*eta; double eta3 = eta*eta2; double eta4 = eta*eta3; switch( NumBasisOrder ) { case 5: // fifth order mu1.set(15, eta4 ); mu1.set(14, xi4 ); mu1.set(13, xi2*eta2 ); mu1.set(12, eta3*xi ); mu1.set(11, xi3*eta ); case 4: // fourth order mu1.set(10, eta3 ); mu1.set(9, xi3 ); mu1.set(8, xi*eta2 ); mu1.set(7, eta*xi2 ); case 3: // third order mu1.set(6, eta2 ); mu1.set(5, xi2 ); mu1.set(4, xi*eta ); case 2: // second order mu1.set(3, eta ); mu1.set(2, xi ); case 1: // first order mu1.set(1, 1.0 ); break; } // Evaluate monomials at locations xy2 xi = xy2.get(1); xi2 = xi*xi; xi3 = xi*xi2; xi4 = xi*xi3; eta = xy2.get(2); eta2 = eta*eta; eta3 = eta*eta2; eta4 = eta*eta3; switch( NumBasisOrder ) { case 5: // fifth order mu2.set(15, eta4 ); mu2.set(14, xi4 ); mu2.set(13, xi2*eta2 ); mu2.set(12, eta3*xi ); mu2.set(11, xi3*eta ); case 4: // fourth order mu2.set(10, eta3 ); mu2.set(9, xi3 ); mu2.set(8, xi*eta2 ); mu2.set(7, eta*xi2 ); case 3: // third order mu2.set(6, eta2 ); mu2.set(5, xi2 ); mu2.set(4, xi*eta ); case 2: // second order mu2.set(3, eta ); mu2.set(2, xi ); case 1: // first order mu2.set(1, 1.0 ); break; } // Finally, convert monomials to Legendre Polys // on the two adjacent triangle for (int k=1; k<=NumBasisComps; k++) { double tmp1 = 0.0; double tmp2 = 0.0; for (int j=1; j<=k; j++) { tmp1 = tmp1 + Mmat[k-1][j-1]*mu1.get(j); tmp2 = tmp2 + Mmat[k-1][j-1]*mu2.get(j); } EdgeData->GL_phi_left->set(i,m,k, tmp1 ); EdgeData->GL_phi_right->set(i,m,k, tmp2 ); } } } }
// Right-hand side for hyperbolic PDE in divergence form // // q_t = -( f(q,x,y,t)_x + g(q,x,y,t)_y ) + Psi(q,x,y,t) // void LaxWendroff_Unst(double dt, const mesh& Mesh, const edge_data_Unst& EdgeData, dTensor3& aux, // SetBndValues modifies ghost cells dTensor3& q, // SetBndValues modifies ghost cells dTensor3& Lstar, dTensor1& smax) { const int NumElems = Mesh.get_NumElems(); const int NumPhysElems = Mesh.get_NumPhysElems(); const int NumEdges = Mesh.get_NumEdges(); const int meqn = q.getsize(2); const int kmax = q.getsize(3); const int maux = aux.getsize(2); const int space_order = dogParams.get_space_order(); dTensor3 EdgeFluxIntegral(NumElems,meqn,kmax); dTensor3 ElemFluxIntegral(NumElems,meqn,kmax); dTensor3 Psi(NumElems,meqn,kmax); // --------------------------------------------------------- // Boundary Conditions SetBndValues_Unst(Mesh, &q, &aux); // --------------------------------------------------------- // --------------------------------------------------------------------- // // Part 0: Compute the Lax-Wendroff "flux" function: // // Here, we include the extra information about time derivatives. // --------------------------------------------------------------------- // dTensor3 F(NumElems, meqn, kmax ); F.setall(0.); dTensor3 G(NumElems, meqn, kmax ); G.setall(0.); L2ProjectLxW_Unst( dogParams.get_time_order(), 1.0, 0.5*dt, dt*dt/6.0, 1, NumElems, space_order, space_order, space_order, space_order, Mesh, &q, &aux, &F, &G, &FluxFunc, &DFluxFunc, &D2FluxFunc ); // --------------------------------------------------------- // Part I: compute source term // --------------------------------------------------------- if ( dogParams.get_source_term()>0 ) { // eprintf("error: have not implemented source term for LxW solver."); printf("Source term has not been implemented for LxW solver. Terminating program."); exit(1); } Lstar.setall(0.); // --------------------------------------------------------- // --------------------------------------------------------- // Part II: compute flux integral on element edges // --------------------------------------------------------- // Loop over all interior edges EdgeFluxIntegral.setall(0.); ElemFluxIntegral.setall(0.); #pragma omp parallel for // Loop over all interior edges for (int i=1; i<=NumEdges; i++) { // Edge coordinates double x1 = Mesh.get_edge(i,1); double y1 = Mesh.get_edge(i,2); double x2 = Mesh.get_edge(i,3); double y2 = Mesh.get_edge(i,4); // Elements on either side of edge int ileft = Mesh.get_eelem(i,1); int iright = Mesh.get_eelem(i,2); double Areal = Mesh.get_area_prim(ileft); double Arear = Mesh.get_area_prim(iright); // Scaled normal to edge dTensor1 nhat(2); nhat.set(1, (y2-y1) ); nhat.set(2, (x1-x2) ); // Variables to store flux integrals along edge dTensor2 Fr_tmp(meqn,dogParams.get_space_order()); dTensor2 Fl_tmp(meqn,dogParams.get_space_order()); // Loop over number of quadrature points along each edge for (int ell=1; ell<=dogParams.get_space_order(); ell++) { dTensor1 Ql(meqn), Qr(meqn); dTensor1 ffl(meqn), ffr(meqn); // << -- NEW PART -- >> dTensor1 Auxl(maux), Auxr(maux); // Riemann data - q for (int m=1; m<=meqn; m++) { Ql.set(m, 0.0 ); Qr.set(m, 0.0 ); // << -- NEW PART, ffl and ffr -- >> // ffl.set(m, 0.0 ); ffr.set(m, 0.0 ); for (int k=1; k<=kmax; k++) { Ql.set(m, Ql.get(m) + EdgeData.phi_left->get(i,ell,k) *q.get(ileft, m,k) ); Qr.set(m, Qr.get(m) + EdgeData.phi_right->get(i,ell,k) *q.get(iright,m,k) ); // << -- NEW PART, ffl and ffr -- >> // // Is this the correct way to use the normal vector? ffl.set(m, ffl.get(m) + EdgeData.phi_left->get (i, ell, k) * ( nhat.get(1)*F.get( ileft, m, k) + nhat.get(2)*G.get( ileft, m, k) ) ); ffr.set(m, ffr.get(m) + EdgeData.phi_right->get(i, ell, k) * ( nhat.get(1)*F.get(iright, m, k) + nhat.get(2)*G.get(iright, m, k) ) ); } } // Riemann data - aux for (int m=1; m<=maux; m++) { Auxl.set(m, 0.0 ); Auxr.set(m, 0.0 ); for (int k=1; k<=kmax; k++) { Auxl.set(m, Auxl.get(m) + EdgeData.phi_left->get(i,ell,k) * aux.get(ileft, m,k) ); Auxr.set(m, Auxr.get(m) + EdgeData.phi_right->get(i,ell,k) * aux.get(iright,m,k) ); } } // Solve Riemann problem dTensor1 xedge(2); double s = EdgeData.xpts1d->get(ell); xedge.set(1, x1 + 0.5*(s+1.0)*(x2-x1) ); xedge.set(2, y1 + 0.5*(s+1.0)*(y2-y1) ); // Solve the Riemann problem for this edge dTensor1 Fl(meqn), Fr(meqn); // Use the time-averaged fluxes to define left and right values for // the Riemann solver. const double smax_edge = RiemannSolveLxW( nhat, xedge, Ql, Qr, Auxl, Auxr, ffl, ffr, Fl, Fr); smax.set(i, Max(smax_edge,smax.get(i)) ); // Construct fluxes for (int m=1; m<=meqn; m++) { Fr_tmp.set(m,ell, Fr.get(m) ); Fl_tmp.set(m,ell, Fl.get(m) ); } } // Add edge integral to line integral around the full element for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { double Fl_sum = 0.0; double Fr_sum = 0.0; for (int ell=1; ell<=dogParams.get_space_order(); ell++) { Fl_sum = Fl_sum + 0.5*EdgeData.wgts1d->get(ell) *EdgeData.phi_left->get(i,ell,k) *Fl_tmp.get(m,ell); Fr_sum = Fr_sum + 0.5*EdgeData.wgts1d->get(ell) *EdgeData.phi_right->get(i,ell,k)*Fr_tmp.get(m,ell); } EdgeFluxIntegral.set(ileft, m,k, EdgeFluxIntegral.get(ileft, m,k) + Fl_sum/Areal ); EdgeFluxIntegral.set(iright,m,k, EdgeFluxIntegral.get(iright,m,k) - Fr_sum/Arear ); } } // --------------------------------------------------------- // --------------------------------------------------------- // Part III: compute intra-element contributions // --------------------------------------------------------- if( dogParams.get_space_order() > 1 ) { L2ProjectGradAddLegendre_Unst(1, NumPhysElems, space_order, Mesh, &F, &G, &ElemFluxIntegral ); } // --------------------------------------------------------- // --------------------------------------------------------- // Part IV: construct Lstar // --------------------------------------------------------- if (dogParams.get_source_term()==0) // Without Source Term { #pragma omp parallel for for (int i=1; i<=NumPhysElems; i++) for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { double tmp = ElemFluxIntegral.get(i,m,k) - EdgeFluxIntegral.get(i,m,k); Lstar.set(i,m,k, tmp ); } } else // With Source Term { #pragma omp parallel for for (int i=1; i<=NumPhysElems; i++) for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { // double tmp = ElemFluxIntegral.get(i,m,k) // - EdgeFluxIntegral.get(i,m,k) // + Psi.get(i,m,k); // Lstar.set(i,m,k, tmp ); printf("Source term has not been implemented for LxW solver. Terminating program."); exit(1); } } // --------------------------------------------------------- // --------------------------------------------------------- // Part V: add extra contributions to Lstar // --------------------------------------------------------- // Call LstarExtra LstarExtra_Unst(Mesh, &q, &aux, &Lstar); // --------------------------------------------------------- // --------------------------------------------------------- // Part VI: artificial viscosity limiter // --------------------------------------------------------- // if (dogParams.get_space_order()>1 && // dogParams.using_viscosity_limiter()) // { ArtificialViscosity(&aux,&q,&Lstar); } // --------------------------------------------------------- }
void ConstructL_Unst( const double t, const dTensor2* vel_vec, const mesh& Mesh, const edge_data_Unst& EdgeData, dTensor3& aux, // SetBndValues_Unst modifies ghost cells dTensor3& q, // SetBndValues_Unst modifies ghost cells dTensor3& Lstar, dTensor1& smax) { const int NumElems = Mesh.get_NumElems(); const int NumPhysElems = Mesh.get_NumPhysElems(); const int NumEdges = Mesh.get_NumEdges(); const int meqn = q.getsize(2); const int kmax = q.getsize(3); const int maux = aux.getsize(2); const int space_order = dogParams.get_space_order(); dTensor3 EdgeFluxIntegral(NumElems,meqn,kmax); dTensor3 ElemFluxIntegral(NumElems,meqn,kmax); dTensor3 Psi(NumElems,meqn,kmax); // --------------------------------------------------------- // Boundary Conditions SetBndValues_Unst(Mesh,&q,&aux); // Positivity limiter void ApplyPosLimiter_Unst(const mesh& Mesh, const dTensor3& aux, dTensor3& q); if( dogParams.using_moment_limiter() ) { ApplyPosLimiter_Unst(Mesh, aux, q); } // --------------------------------------------------------- // --------------------------------------------------------- // Part I: compute flux integral on element edges // --------------------------------------------------------- // Loop over all interior edges and solve Riemann problems // dTensor1 nvec(2); // Loop over all interior edges EdgeFluxIntegral.setall(0.); ElemFluxIntegral.setall(0.); // Loop over all interior edges #pragma omp parallel for for (int i=1; i<=NumEdges; i++) { // Edge coordinates double x1 = Mesh.get_edge(i,1); double y1 = Mesh.get_edge(i,2); double x2 = Mesh.get_edge(i,3); double y2 = Mesh.get_edge(i,4); // Elements on either side of edge int ileft = Mesh.get_eelem(i,1); int iright = Mesh.get_eelem(i,2); double Areal = Mesh.get_area_prim(ileft); double Arear = Mesh.get_area_prim(iright); // Scaled normal to edge dTensor1 nhat(2); nhat.set(1, (y2-y1) ); nhat.set(2, (x1-x2) ); // Variables to store flux integrals along edge dTensor2 Fr_tmp(meqn,dogParams.get_space_order()); dTensor2 Fl_tmp(meqn,dogParams.get_space_order()); // Loop over number of quadrature points along each edge for (int ell=1; ell<=dogParams.get_space_order(); ell++) { dTensor1 Ql(meqn),Qr(meqn); dTensor1 Auxl(maux),Auxr(maux); // Riemann data - q for (int m=1; m<=meqn; m++) { Ql.set(m, 0.0 ); Qr.set(m, 0.0 ); for (int k=1; k<=kmax; k++) { Ql.set(m, Ql.get(m) + EdgeData.phi_left->get(i,ell,k) *q.get(ileft, m,k) ); Qr.set(m, Qr.get(m) + EdgeData.phi_right->get(i,ell,k) *q.get(iright,m,k) ); } } // Riemann data - aux for (int m=1; m<=maux; m++) { Auxl.set(m, 0.0 ); Auxr.set(m, 0.0 ); for (int k=1; k<=kmax; k++) { Auxl.set(m, Auxl.get(m) + EdgeData.phi_left->get(i,ell,k) *aux.get(ileft, m,k) ); Auxr.set(m, Auxr.get(m) + EdgeData.phi_right->get(i,ell,k) *aux.get(iright,m,k) ); } } // Solve Riemann problem dTensor1 xedge(2); double s = EdgeData.xpts1d->get(ell); xedge.set(1, x1 + 0.5*(s+1.0)*(x2-x1) ); xedge.set(2, y1 + 0.5*(s+1.0)*(y2-y1) ); dTensor1 Fl(meqn),Fr(meqn); const double smax_edge = RiemannSolve(vel_vec, nhat, xedge, Ql, Qr, Auxl, Auxr, Fl, Fr); smax.set(i, Max(smax_edge,smax.get(i)) ); // Construct fluxes for (int m=1; m<=meqn; m++) { Fr_tmp.set(m,ell, Fr.get(m) ); Fl_tmp.set(m,ell, Fl.get(m) ); } } // Add edge integral to line integral around the full element for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { double Fl_sum = 0.0; double Fr_sum = 0.0; for (int ell=1; ell<=dogParams.get_space_order(); ell++) { Fl_sum = Fl_sum + 0.5*EdgeData.wgts1d->get(ell) *EdgeData.phi_left->get(i,ell,k) *Fl_tmp.get(m,ell); Fr_sum = Fr_sum + 0.5*EdgeData.wgts1d->get(ell) *EdgeData.phi_right->get(i,ell,k)*Fr_tmp.get(m,ell); } EdgeFluxIntegral.set(ileft, m,k, EdgeFluxIntegral.get(ileft, m,k) + Fl_sum/Areal ); EdgeFluxIntegral.set(iright,m,k, EdgeFluxIntegral.get(iright,m,k) - Fr_sum/Arear ); } } // --------------------------------------------------------- // --------------------------------------------------------- // Part II: compute intra-element contributions // --------------------------------------------------------- L2ProjectGrad_Unst(vel_vec, 1,NumPhysElems, space_order,space_order,space_order,space_order, Mesh,&q,&aux,&ElemFluxIntegral,&FluxFunc); // --------------------------------------------------------- // --------------------------------------------------------- // Part III: compute source term // --------------------------------------------------------- if ( dogParams.get_source_term()>0 ) { // Set source term on computational grid // Set values and apply L2-projection L2Project_Unst(t, vel_vec, 1,NumPhysElems, space_order,space_order,space_order,space_order, Mesh,&q,&aux,&Psi,&SourceTermFunc); } // --------------------------------------------------------- // --------------------------------------------------------- // Part IV: construct Lstar // --------------------------------------------------------- if (dogParams.get_source_term()==0) // Without Source Term { #pragma omp parallel for for (int i=1; i<=NumPhysElems; i++) for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { double tmp = ElemFluxIntegral.get(i,m,k) - EdgeFluxIntegral.get(i,m,k); Lstar.set(i,m,k, tmp ); } } else // With Source Term { #pragma omp parallel for for (int i=1; i<=NumPhysElems; i++) for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { double tmp = ElemFluxIntegral.get(i,m,k) - EdgeFluxIntegral.get(i,m,k) + Psi.get(i,m,k); Lstar.set(i,m,k, tmp ); } } // --------------------------------------------------------- // --------------------------------------------------------- // Part V: add extra contributions to Lstar // --------------------------------------------------------- // Call LstarExtra LstarExtra_Unst(Mesh,&q,&aux,&Lstar); // --------------------------------------------------------- }