// This is a user-supplied routine that sets the the boundary conditions // void SetBndValues_Unst(const mesh& Mesh, dTensor3* q, dTensor3* aux) { int meqn = q->getsize(2); int kmax = q->getsize(3); int maux = aux->getsize(2); int NumElems = Mesh.get_NumElems(); int NumPhysElems = Mesh.get_NumPhysElems(); int NumGhostElems = Mesh.get_NumGhostElems(); int NumNodes = Mesh.get_NumNodes(); int NumPhysNodes = Mesh.get_NumPhysNodes(); int NumEdges = Mesh.get_NumEdges(); // ---------------------------------------- // Loop over each ghost cell element and // place the correct information into // these elements // ---------------------------------------- for (int i=1; i<=NumGhostElems; i++) { int j = Mesh.get_ghost_link(i); for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { q->set(i+NumPhysElems,m,k, 0.0 ); } } }
double GetCFL_Unst(double dt, const mesh& Mesh, const dTensor3& aux, const dTensor1& smax) { double cfl=-100.0; int NumPhysElems = Mesh.get_NumPhysElems(); int NumEdges = Mesh.get_NumEdges(); for (int i=1; i<=NumPhysElems; i++) { double Area = Mesh.get_area_prim(i); int edge1 = Mesh.get_tedge(i,1); int edge2 = Mesh.get_tedge(i,2); int edge3 = Mesh.get_tedge(i,3); double tmp = Max(Max(smax.get(edge1),smax.get(edge2)), smax.get(edge3)); cfl = Max(0.5*dt*tmp/Area, cfl); } return cfl; }
// This is a user-supplied routine that sets the the boundary conditions // // The default routine for the 4D Vlasov code is to apply zero boundary // conditions in configuration space. This routine is identical to the one // found in the unst branch of the 2D DoGPack code, and *should* be setting // periodic boundary conditions. // void SetBndValues_Unst(const mesh& Mesh, dTensor3* q, dTensor3* aux) { // problem information (If this were to be pulled from DogParams, these // numbers would NOT be correct! The reason is that each quadrature point // was actually saved as a separate "equation") const int meqn = q->getsize(2); const int kmax = q->getsize(3); const int maux = aux->getsize(2); // Mesh information const int NumElems = Mesh.get_NumElems(); const int NumPhysElems = Mesh.get_NumPhysElems(); const int NumGhostElems = Mesh.get_NumGhostElems(); const int NumNodes = Mesh.get_NumNodes(); const int NumPhysNodes = Mesh.get_NumPhysNodes(); const int NumEdges = Mesh.get_NumEdges(); // ---------------------------------------- // Loop over each ghost cell element and // place the correct information into // these elements // ---------------------------------------- for (int i=1; i<=NumGhostElems; i++) { int j = Mesh.get_ghost_link(i); for (int m=1; m<=meqn; m++) for (int k=1; k<=kmax; k++) { q->set(i+NumPhysElems, m, k, q->get(j, m, k) ); } } }
// // Output basic mesh information to screen // void ScreenOutput(const mesh& Mesh) { // Compute mesh quality parameters double totalarea = Mesh.get_area_prim(1); double maxarea = Mesh.get_area_prim(1); double minarea = Mesh.get_area_prim(1); for (int i=2; i<=Mesh.get_NumPhysElems(); i++) { double tmp = Mesh.get_area_prim(i); totalarea = totalarea + tmp; if (tmp < minarea) { minarea = tmp; } if (tmp > maxarea) { maxarea = tmp; } } double minAngle = 180.0; for (int i=1; i<=Mesh.get_NumPhysElems(); i++) { const int i1 = Mesh.get_tnode(i,1); const int i2 = Mesh.get_tnode(i,2); const int i3 = Mesh.get_tnode(i,3); point v12, v23, v31; v12.x = Mesh.get_node(i2,1) - Mesh.get_node(i1,1); v12.y = Mesh.get_node(i2,2) - Mesh.get_node(i1,2); v23.x = Mesh.get_node(i3,1) - Mesh.get_node(i2,1); v23.y = Mesh.get_node(i3,2) - Mesh.get_node(i2,2); v31.x = Mesh.get_node(i1,1) - Mesh.get_node(i3,1); v31.y = Mesh.get_node(i1,2) - Mesh.get_node(i3,2); double angle1 = acos((v12.x*-v31.x+v12.y*-v31.y) /(sqrt(v12.x*v12.x+v12.y*v12.y)*sqrt(v31.x*v31.x+v31.y*v31.y))); double angle2 = acos((v23.x*-v12.x+v23.y*-v12.y) /(sqrt(v23.x*v23.x+v23.y*v23.y)*sqrt(v12.x*v12.x+v12.y*v12.y))); double angle3 = acos((v31.x*-v23.x+v31.y*-v23.y) /(sqrt(v31.x*v31.x+v31.y*v31.y)*sqrt(v23.x*v23.x+v23.y*v23.y))); if ((angle1*180/pi) < minAngle) { minAngle = angle1*180/pi; } if ((angle2*180/pi) < minAngle) { minAngle = angle2*180/pi; } if ((angle3*180/pi) < minAngle) { minAngle = angle3*180/pi; } } // Output summary of results to screen printf("\n"); printf(" SUMMARY OF RESULTS:\n"); printf(" -------------------\n"); printf(" Number of Elements: %8i\n",Mesh.get_NumElems()); printf(" Number of Physical Elements: %8i\n",Mesh.get_NumPhysElems()); printf(" Number of Ghost Elements: %8i\n",Mesh.get_NumGhostElems()); printf(" Number of Nodes: %8i\n",Mesh.get_NumNodes()); printf(" Number of Physical Nodes: %8i\n",Mesh.get_NumPhysNodes()); printf(" Number of Boundary Nodes: %8i\n",Mesh.get_NumBndNodes()); printf(" Number of Edges: %8i\n",Mesh.get_NumEdges()); printf(" Number of Boundary Edges: %8i\n",Mesh.get_NumBndEdges()); printf("\n"); printf(" Total Area Covered: %24.16e\n",totalarea); printf(" Area Ratio: small/large: %24.16e\n",minarea/maxarea); printf(" Angle Ratio: minAngle/60: %24.16e\n",minAngle/60.0); printf("\n"); }
// 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); } // --------------------------------------------------------- }
// This file should be identical to DogSolveRK_Unst, with the exception that all // output printing statements are silenced. // // Advance the solution qold to qnew over time interval tstart to tend. // // All local information is allocated within this function. The only part // that gets shared are time values passed through dogStateUnst2. This class // should be modified to accept the state variable, q and aux in place of only // containing time information as is currently the case. (-DS) double DogSolveRK_Unst_Quiet( const dTensor2* vel_vec, const mesh& Mesh, const edge_data_Unst& EdgeData, dTensor3& aux, dTensor3& qold, dTensor3& qnew, const double tstart, const double tend, DogStateUnst2& dogStateUnst2) { const int mx = qnew.getsize(1); const int meqn = qnew.getsize(2); const int kmax = qnew.getsize(3); const int maux = aux.getsize(2); const double* cflv = dogParams.get_cflv(); const int nv = dogParams.get_nv(); RKinfo rk; SetRKinfo(dogParams.get_time_order(),rk); // define local variables int n_step = 0; double t = tstart; double dt = dogStateUnst2.get_initial_dt(); const double CFL_max = cflv[1]; const double CFL_target = cflv[2]; double cfl = 0.0; double dtmin = dt; double dtmax = dt; const int NumElems = Mesh.get_NumElems(); // Number of total elements in mesh const int NumNodes = Mesh.get_NumNodes(); // Number of nodes in mesh const int NumEdges = Mesh.get_NumEdges(); // Number of edges in mesh dTensor3 qstar(NumElems,meqn,kmax); dTensor3 q1(NumElems,meqn,kmax); dTensor3 q2(NumElems,meqn,kmax); dTensor3 auxstar(NumElems,maux,kmax); dTensor3 Lstar(NumElems,meqn,kmax); dTensor3 Lold(NumElems,meqn,kmax); dTensor3 auxold(NumElems,maux,kmax); dTensor1 smax(NumEdges); void L2Project_Unst( const dTensor2* vel_vec, const int istart, const int iend, const int QuadOrder, const int BasisOrder_qin, const int BasisOrder_auxin, const int BasisOrder_fout, const mesh& Mesh, const dTensor3* qin, const dTensor3* auxin, dTensor3* fout, void (*Func)(const dTensor2* vel_vec, const dTensor2&,const dTensor2&, const dTensor2&,dTensor2&)); // JUNK here: void AuxFuncWrapper( const dTensor2* vel_vec, const dTensor2& xpts, const dTensor2& NOT_USED_1, const dTensor2& NOT_USED_2, dTensor2& auxvals); const int space_order = dogParams.get_space_order(); if( maux > 0 ) { printf("WARNING: maux = %d should be zero for Vlasov routines.", maux); printf(" Modify parameters.ini to remove this warning\n" ); L2Project_Unst(vel_vec,1,NumElems, space_order,space_order,space_order,space_order, Mesh,&qnew,&aux,&aux,&AuxFuncWrapper); } // Set initialize qstar and auxstar values qstar.copyfrom(qold); auxstar.copyfrom(aux); // Runge-Kutta time stepping while (t<tend) { // initialize time step int m_accept = 0; n_step = n_step + 1; // check if max number of time steps exceeded if( n_step > nv ) { eprintf(" Error in DogSolveRK_Unst.cpp: " " Exceeded allowed # of time steps \n" " n_step = %d\n" " nv = %d\n\n", n_step, nv); } // copy qnew into qold qold.copyfrom(qnew); auxold.copyfrom(aux); // keep trying until we get a dt that does not violate CFL condition while (m_accept==0) { // set current time double told = t; if (told+dt > tend) { dt = tend - told; } t = told + dt; // TODO - this needs to be performed at the 'local' level dogStateUnst2.set_time ( told ); dogStateUnst2.set_dt ( dt ); // Set initial maximum wave speed to zero smax.setall(0.); // Take a full time step of size dt switch ( dogParams.get_time_order() ) { case 1: // First order in time (1-stage) // ----------------------------------------------- // Stage #1 (the only one in this case) rk.mstage = 1; BeforeStep_Unst(dt,Mesh,aux,qnew); ConstructL_Unst(told, vel_vec,Mesh,EdgeData,aux,qnew,Lstar,smax); UpdateSoln_Unst(rk.alpha1->get(rk.mstage),rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage),dt,Mesh,aux, qnew, Lstar, qnew); AfterStep_Unst(dt,Mesh,aux,qnew); // ----------------------------------------------- break; case 2: // Second order in time (2-stages) // ----------------------------------------------- // Stage #1 rk.mstage = 1; dogStateUnst2.set_time(told); BeforeStep_Unst(dt,Mesh,aux,qnew); ConstructL_Unst(told,vel_vec,Mesh,EdgeData,aux,qnew,Lstar,smax); UpdateSoln_Unst( rk.alpha1->get(rk.mstage),rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage), dt, Mesh, aux, qnew, Lstar, qstar); AfterStep_Unst(dt, Mesh, auxstar, qstar); // ------------------------------------------------ // Stage #2 rk.mstage = 2; dogStateUnst2.set_time(told+dt); BeforeStep_Unst(dt, Mesh, auxstar, qstar); ConstructL_Unst(told+1.0*dt, vel_vec, Mesh, EdgeData, aux, qstar, Lstar, smax); UpdateSoln_Unst(rk.alpha1->get(rk.mstage), rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage), dt, Mesh, auxstar, qstar, Lstar, qnew); AfterStep_Unst(dt, Mesh, aux, qnew); // ------------------------------------------------ break; case 3: // Third order in time (3-stages) // qnew = alpha1 * qstar + alpha2 * qnew + beta * dt * L( qstar ) // alpha1 = 1.0 // alpha2 = 0.0 // beta = 1.0 // ------------------------------------------------ // Stage #1 rk.mstage = 1; dogStateUnst2.set_time(told); BeforeStep_Unst(dt,Mesh,aux,qnew); ConstructL_Unst(told, vel_vec,Mesh,EdgeData,aux,qnew,Lstar,smax); Lold.copyfrom(Lstar); UpdateSoln_Unst(rk.alpha1->get(rk.mstage),rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage),dt,Mesh,aux,qnew,Lstar,qstar); AfterStep_Unst(dt,Mesh,aux,qstar); // ------------------------------------------------- // alpha1 = 0.75 // alpha2 = 0.25 // beta = 0.25 // Stage #2 rk.mstage = 2; dogStateUnst2.set_time(told+0.5*dt); BeforeStep_Unst(dt,Mesh,aux,qstar); ConstructL_Unst(told+dt, vel_vec,Mesh,EdgeData,aux,qstar,Lstar,smax); UpdateSoln_Unst(rk.alpha1->get(rk.mstage),rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage),dt,Mesh,aux,qnew,Lstar,qstar); AfterStep_Unst(dt,Mesh,aux,qstar); // -------------------------------------------------- // alpha1 = 2/3 // alpha2 = 1/3 // beta = 2/3 // Stage #3 rk.mstage = 3; dogStateUnst2.set_time(told+dt); BeforeStep_Unst(dt,Mesh,auxstar,qstar); ConstructL_Unst(told+0.5*dt,vel_vec,Mesh,EdgeData,auxstar,qstar,Lstar,smax); UpdateSoln_Unst(rk.alpha1->get(rk.mstage),rk.alpha2->get(rk.mstage), rk.beta->get(rk.mstage),dt,Mesh,aux,qstar,Lstar,qnew); AfterStep_Unst(dt,Mesh,aux,qnew); // -------------------------------------------------- break; default: unsupported_value_error(dogParams.get_time_order()); } // compute cfl number cfl = GetCFL_Unst(dt,Mesh,aux,smax); // output time step information // if (dogParams.get_verbosity()>0) // { // printf(" In DogSolveRK_Quiet: DogSolve2D ... Step %5d" // " CFL =%6.3f" // " dt =%11.3e" // " t =%11.3e\n", // n_step, cfl, dt, t); // } // choose new time step if (cfl>0.0) { dt = Min(dogParams.get_max_dt(), dt*CFL_target/cfl); dtmin = Min(dt,dtmin); dtmax = Max(dt,dtmax); } else { dt = dogParams.get_max_dt(); } // see whether to accept or reject this step if (cfl<=CFL_max) // accept { m_accept = 1; dogStateUnst2.set_time(t); // do any extra work // AfterFullTimeStep_Unst(dogStateUnst2.get_dt(),Mesh, // auxold,qold,Lold,aux,qnew); } else //reject { t = told; dogStateUnst2.set_time(told); // if( dogParams.get_verbosity() > 0 ) // { // printf("DogSolve2D rejecting step..." // "CFL number too large\n"); // } // copy qold into qnew qnew.copyfrom(qold); aux.copyfrom(auxold); // after reject function // AfterReject_Unst(Mesh,dt,aux,qnew); } } } // printf(" Finished! t = %2.3e and nsteps = %d\n", t, n_step ); // set initial time step for next call to DogSolveRK dogStateUnst2.set_initial_dt(dt); void DeleteRKInfo(RKinfo& rk); DeleteRKInfo(rk); return cfl; }
// This is the positivity preserving limiter proposed in // "Maximum-Principle-Satisfying and Positivity-Preserving // High Order Discontinuous Galerkin Schemes // for Conservation Laws on Triangular Meshes", Zhang, Xia and Shu // J. Sci. Comput. (2012). // // THIS METHOD ASSUMES THAT EVERY COMPONENT OF CONSERVED VARIABLES SHOULD STAY // POSITIVE. // // In order to implement this for a different scheme, one should rewrite, or // redefine what components should remain positiive. This will require // reworking the control flow logic for how time step lengths are chosen. void ApplyPosLimiter_Unst(const mesh& Mesh, const dTensor3& aux, dTensor3& q) { 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(); // Do nothing in the case of piecewise constants if( space_order == 1 ) { return; } // ------------------------------------------------ // // number of points where we want to check solution // // ------------------------------------------------ // const int space_order_sq = space_order*space_order; const int mpts_vec[] = {0, 3*space_order_sq, 18, 3*space_order_sq, 3*space_order_sq }; // TODO - FILL IN 2ND-ORDER CASE const int mpoints = mpts_vec[space_order-1]; // ---------------------------------------------------------- // // sample basis at all points where we want to check solution // // ---------------------------------------------------------- // dTensor2 spts(mpoints, 2); void SetPositivePoints_Unst(const int& space_order, dTensor2& spts); SetPositivePoints_Unst(space_order, spts); void SamplePhiAtPositivePoints_Unst(const int& space_order, const dTensor2& spts, dTensor2& phi); dTensor2 phi(mpoints, kmax); SamplePhiAtPositivePoints_Unst(space_order, spts, phi); // -------------------------------------------------------------- // // q_limited = Q1 + \theta ( q(xi,eta) - Q1 ) // // where theta = min(1, |Q1| / |Q1-m|; m = min_{i} q(xi_i, eta_i) // // -------------------------------------------------------------- // #pragma omp parallel for for(int i=1; i <= NumPhysElems; i++) for(int me=1; me <= meqn; me++) { double m = 0.0; for(int mp=1; mp <= mpoints; mp++) { // evaluate q at spts(mp) // double qnow = 0.0; for( int k=1; k <= kmax; k++ ) { qnow += q.get(i,me,k) * phi.get(mp,k); } m = Min(m, qnow); } double theta = 0.0; double Q1 = q.get(i,me,1); assert_ge( Q1, -1e-13 ); if( fabs( Q1 - m ) < 1.0e-14 ){ theta = 1.0; } else{ theta = Min( 1.0, fabs( Q1 / (Q1 - m) ) ); } // limit q // for( int k=2; k <= kmax; k++ ) { q.set(i,me,k, q.get(i,me,k) * theta ); } } }
// 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 ); } } } }
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); // --------------------------------------------------------- }