void EulerADCFE::DoOdeProjection( const Array<OneD, const Array<OneD, NekDouble> > &inarray, Array<OneD, Array<OneD, NekDouble> > &outarray, const NekDouble time) { int i; int nvariables = inarray.num_elements(); switch(m_projectionType) { case MultiRegions::eDiscontinuous: { // Just copy over array int npoints = GetNpoints(); for(i = 0; i < nvariables; ++i) { Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1); } SetBoundaryConditions(outarray, time); break; } case MultiRegions::eGalerkin: case MultiRegions::eMixed_CG_Discontinuous: { ASSERTL0(false, "No Continuous Galerkin for Euler equations"); break; } default: ASSERTL0(false, "Unknown projection scheme"); break; } }
void EigenValuesAdvection::v_GetFluxVector(const int i, Array<OneD, Array<OneD, NekDouble> > &physfield, Array<OneD, Array<OneD, NekDouble> > &flux) { ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match"); for(int j = 0; j < flux.num_elements(); ++j) { Vmath::Vmul(GetNpoints(),physfield[i],1, m_velocity[j],1,flux[j],1); } }
void IncNavierStokes::v_GetFluxVector(const int i, Array<OneD, Array<OneD, NekDouble> > &physfield, Array<OneD, Array<OneD, NekDouble> > &flux) { ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match"); for(int j = 0; j < flux.num_elements(); ++j) { Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1); } }
/** * @brief Compute the projection for the linear advection equation. * * @param inarray Given fields. * @param outarray Calculated solution. * @param time Time. */ void UnsteadyAdvection::DoOdeProjection( const Array<OneD, const Array<OneD, NekDouble> >&inarray, Array<OneD, Array<OneD, NekDouble> >&outarray, const NekDouble time) { // Counter variable int i; // Number of fields (variables of the problem) int nVariables = inarray.num_elements(); // Set the boundary conditions SetBoundaryConditions(time); // Switch on the projection type (Discontinuous or Continuous) switch(m_projectionType) { // Discontinuous projection case MultiRegions::eDiscontinuous: { // Number of quadrature points int nQuadraturePts = GetNpoints(); // Just copy over array for(i = 0; i < nVariables; ++i) { Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1); } break; } // Continuous projection case MultiRegions::eGalerkin: case MultiRegions::eMixed_CG_Discontinuous: { Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(),0.0); for(i = 0; i < nVariables; ++i) { m_fields[i]->FwdTrans(inarray[i], coeffs); m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]); } break; } default: ASSERTL0(false,"Unknown projection scheme"); break; } }
/** * @brief Compute the right-hand side for the linear advection equation. * * @param inarray Given fields. * @param outarray Calculated solution. * @param time Time. */ void UnsteadyAdvection::DoOdeRhs( const Array<OneD, const Array<OneD, NekDouble> >&inarray, Array<OneD, Array<OneD, NekDouble> >&outarray, const NekDouble time) { // Counter variable int i; // Number of fields (variables of the problem) int nVariables = inarray.num_elements(); // Number of solution points int nSolutionPts = GetNpoints(); // RHS computation using the new advection base class m_advection->Advect(nVariables, m_fields, m_velocity, inarray, outarray); // Negate the RHS for (i = 0; i < nVariables; ++i) { Vmath::Neg(nSolutionPts, outarray[i], 1); } }
void EigenValuesAdvection::v_DoSolve() { int nvariables = 1; int i,dofs = GetNcoeffs(); //bool UseContCoeffs = false; Array<OneD, Array<OneD, NekDouble> > inarray(nvariables); Array<OneD, Array<OneD, NekDouble> > tmp(nvariables); Array<OneD, Array<OneD, NekDouble> > outarray(nvariables); Array<OneD, Array<OneD, NekDouble> > WeakAdv(nvariables); int npoints = GetNpoints(); int ncoeffs = GetNcoeffs(); switch (m_projectionType) { case MultiRegions::eDiscontinuous: { dofs = ncoeffs; break; } case MultiRegions::eGalerkin: case MultiRegions::eMixed_CG_Discontinuous: { //dofs = GetContNcoeffs(); //UseContCoeffs = true; break; } } cout << endl; cout << "Num Phys Points = " << npoints << endl; // phisical points cout << "Num Coeffs = " << ncoeffs << endl; // cout << "Num Cont Coeffs = " << dofs << endl; inarray[0] = Array<OneD, NekDouble>(npoints,0.0); outarray[0] = Array<OneD, NekDouble>(npoints,0.0); tmp[0] = Array<OneD, NekDouble>(npoints,0.0); WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs,0.0); Array<OneD, NekDouble> MATRIX(npoints*npoints,0.0); for (int j = 0; j < npoints; j++) { inarray[0][j] = 1.0; /// Feeding the weak Advection oprator with a vector (inarray) /// Looping on inarray and changing the position of the only non-zero entry /// we simulate the multiplication by the identity matrix. /// The results stored in outarray is one of the columns of the weak advection oprators /// which are then stored in MATRIX for the futher eigenvalues calculation. switch (m_projectionType) { case MultiRegions::eDiscontinuous: { WeakDGAdvection(inarray, WeakAdv,true,true,1); m_fields[0]->MultiplyByElmtInvMass(WeakAdv[0],WeakAdv[0]); m_fields[0]->BwdTrans(WeakAdv[0],outarray[0]); Vmath::Neg(npoints,outarray[0],1); break; } case MultiRegions::eGalerkin: case MultiRegions::eMixed_CG_Discontinuous: { // Calculate -V\cdot Grad(u); for(i = 0; i < nvariables; ++i) { //Projection m_fields[i]->FwdTrans(inarray[i],WeakAdv[i]); m_fields[i]->BwdTrans_IterPerExp(WeakAdv[i],tmp[i]); //Advection operator AdvectionNonConservativeForm(m_velocity,tmp[i],outarray[i]); Vmath::Neg(npoints,outarray[i],1); //m_fields[i]->MultiplyByInvMassMatrix(WeakAdv[i],WeakAdv[i]); //Projection m_fields[i]->FwdTrans(outarray[i],WeakAdv[i]); m_fields[i]->BwdTrans_IterPerExp(WeakAdv[i],outarray[i]); } break; } } /// The result is stored in outarray (is the j-th columns of the weak advection operator). /// We now store it in MATRIX(j) Vmath::Vcopy(npoints,&(outarray[0][0]),1,&(MATRIX[j]),npoints); /// Set the j-th entry of inarray back to zero inarray[0][j] = 0.0; } //////////////////////////////////////////////////////////////////////////////// /// Calulating the eigenvalues of the weak advection operator stored in (MATRIX) /// using Lapack routines char jobvl = 'N'; char jobvr = 'N'; int info = 0, lwork = 3*npoints; NekDouble dum; Array<OneD, NekDouble> EIG_R(npoints); Array<OneD, NekDouble> EIG_I(npoints); Array<OneD, NekDouble> work(lwork); Lapack::Dgeev(jobvl,jobvr,npoints,MATRIX.get(),npoints,EIG_R.get(),EIG_I.get(),&dum,1,&dum,1,&work[0],lwork,info); //////////////////////////////////////////////////////// //Print Matrix FILE *mFile; mFile = fopen ("WeakAdvMatrix.txt","w"); for(int j = 0; j<npoints; j++) { for(int k = 0; k<npoints; k++) { fprintf(mFile,"%e ",MATRIX[j*npoints+k]); } fprintf(mFile,"\n"); } fclose (mFile); //////////////////////////////////////////////////////// //Output of the EigenValues FILE *pFile; pFile = fopen ("Eigenvalues.txt","w"); for(int j = 0; j<npoints; j++) { fprintf(pFile,"%e %e\n",EIG_R[j],EIG_I[j]); } fclose (pFile); cout << "\nEigenvalues : " << endl; for(int j = 0; j<npoints; j++) { cout << EIG_R[j] << "\t" << EIG_I[j] << endl; } cout << endl; }
void EulerADCFE::DoOdeRhs( const Array<OneD, const Array<OneD, NekDouble> > &inarray, Array<OneD, Array<OneD, NekDouble> > &outarray, const NekDouble time) { int i; int nvariables = inarray.num_elements(); int npoints = GetNpoints(); Array<OneD, Array<OneD, NekDouble> > advVel; Array<OneD, Array<OneD, NekDouble> > outarrayAdv(nvariables); Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nvariables); for (i = 0; i < nvariables; ++i) { outarrayAdv[i] = Array<OneD, NekDouble>(npoints, 0.0); outarrayDiff[i] = Array<OneD, NekDouble>(npoints, 0.0); } m_advection->Advect(nvariables, m_fields, advVel, inarray, outarrayAdv, m_time); for (i = 0; i < nvariables; ++i) { Vmath::Neg(npoints, outarrayAdv[i], 1); } m_diffusion->Diffuse(nvariables, m_fields, inarray, outarrayDiff); if (m_shockCaptureType == "NonSmooth") { for (i = 0; i < nvariables; ++i) { Vmath::Vadd(npoints, outarrayAdv[i], 1, outarrayDiff[i], 1, outarray[i], 1); } } if(m_shockCaptureType == "Smooth") { const Array<OneD, int> ExpOrder = GetNumExpModesPerExp(); NekDouble pOrder = Vmath::Vmax(ExpOrder.num_elements(), ExpOrder, 1); Array <OneD, NekDouble > a_vel (npoints, 0.0); Array <OneD, NekDouble > u_abs (npoints, 0.0); Array <OneD, NekDouble > pres (npoints, 0.0); Array <OneD, NekDouble > wave_sp(npoints, 0.0); GetPressure(inarray, pres); GetSoundSpeed(inarray, pres, a_vel); GetAbsoluteVelocity(inarray, u_abs); Vmath::Vadd(npoints, a_vel, 1, u_abs, 1, wave_sp, 1); NekDouble max_wave_sp = Vmath::Vmax(npoints, wave_sp, 1); Vmath::Smul(npoints, m_C2, outarrayDiff[nvariables-1], 1, outarrayDiff[nvariables-1], 1); Vmath::Smul(npoints, max_wave_sp, outarrayDiff[nvariables-1], 1, outarrayDiff[nvariables-1], 1); Vmath::Smul(npoints, pOrder, outarrayDiff[nvariables-1], 1, outarrayDiff[nvariables-1], 1); for (i = 0; i < nvariables; ++i) { Vmath::Vadd(npoints, outarrayAdv[i], 1, outarrayDiff[i], 1, outarray[i], 1); } Array<OneD, Array<OneD, NekDouble> > outarrayForcing(nvariables); for (i = 0; i < nvariables; ++i) { outarrayForcing[i] = Array<OneD, NekDouble>(npoints, 0.0); } GetForcingTerm(inarray, outarrayForcing); for (i = 0; i < nvariables; ++i) { // Add Forcing Term Vmath::Vadd(npoints, outarray[i], 1, outarrayForcing[i], 1, outarray[i], 1); } } // Add sponge layer if defined in the session file std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x; for (x = m_forcing.begin(); x != m_forcing.end(); ++x) { (*x)->Apply(m_fields, inarray, outarray, time); } }