43 namespace LibUtilities
80 m_solVector(m_scheme->GetNsteps()),
81 m_t(m_scheme->GetNsteps())
87 int nvar = y.num_elements();
88 int npoints = y[0].num_elements();
89 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
91 for(
int i = 1; i < nsteps; i++)
93 m_solVector[i] = Array<OneD, Array<OneD, NekDouble> >(nvar);
94 for(
int j = 0; j < nvar; j++)
96 m_solVector[i][j] = Array<OneD,NekDouble>(npoints,0.0);
98 if(i < nMultiStepVals)
100 m_t[i] = time - i*timestep*timeLevels[i];
111 const Array<OneD, NekDouble>& t):
116 ASSERTL1(y.num_elements()==
m_scheme->GetNsteps(),
"Amount of Entries does not match number of (multi-) steps");
121 unsigned int npoints):
123 m_solVector(m_scheme->GetNsteps()),
124 m_t(m_scheme->GetNsteps())
126 for(
int i = 0; i <
m_scheme->GetNsteps(); i++)
128 m_solVector[i] = Array<OneD, Array<OneD, NekDouble> >(nvar);
129 for(
int j = 0; j < nvar; j++)
131 m_solVector[i][j] = Array<OneD,NekDouble>(npoints);
138 m_solVector(m_scheme->GetNsteps()),
139 m_t(m_scheme->GetNsteps())
162 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
163 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
182 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
183 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
190 m_B[0][0][0] = 3.0/2.0;
203 m_timeLevelOffset[1] = 1;
211 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
212 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
222 m_U[0][1] = 23.0/12.0;
223 m_U[0][2] = -4.0/3.0;
224 m_U[0][3] = 5.0/12.0;
227 m_V[0][1] = 23.0/12.0;
228 m_V[0][2] = -4.0/3.0;
229 m_V[0][3] = 5.0/12.0;
238 m_timeLevelOffset[1] = 1;
239 m_timeLevelOffset[2] = 2;
240 m_timeLevelOffset[3] = 3;
250 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
251 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
270 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
271 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
290 m_timeLevelOffset[1] = 0;
299 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
300 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
309 m_B[0][0][0] = 2*third;
312 m_U[0][3] = -2*third;
317 m_V[0][3] = -2*third;
326 m_timeLevelOffset[1] = 1;
327 m_timeLevelOffset[2] = 0;
328 m_timeLevelOffset[3] = 1;
337 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
338 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
347 m_B[0][0][0] = 6*eleventh;
349 m_U[0][1] = -9*eleventh;
350 m_U[0][2] = 2*eleventh;
351 m_U[0][4] = -18*eleventh;
352 m_U[0][5] = 6*eleventh;
354 m_V[0][0] = 18*eleventh;
355 m_V[0][1] = -9*eleventh;
356 m_V[0][2] = 2*eleventh;
357 m_V[0][3] = 18*eleventh;
358 m_V[0][4] = -18*eleventh;
359 m_V[0][5] = 6*eleventh;
370 m_timeLevelOffset[1] = 1;
371 m_timeLevelOffset[2] = 2;
372 m_timeLevelOffset[3] = 0;
373 m_timeLevelOffset[4] = 1;
374 m_timeLevelOffset[5] = 2;
382 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
383 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
405 m_timeLevelOffset[1] = 0;
414 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
415 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
422 m_B[0][0][0] = 2*third;
437 m_timeLevelOffset[1] = 1;
445 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
446 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
468 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
469 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
480 m_B[0][0][0] = 1.0/6.0;
481 m_B[0][0][1] = 1.0/3.0;
482 m_B[0][0][2] = 1.0/3.0;
483 m_B[0][0][3] = 1.0/6.0;
497 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
498 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
522 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
523 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
544 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
545 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
554 m_A[0][0][0] = lambda;
555 m_A[0][1][0] = 1.0 - lambda;
556 m_A[0][1][1] = lambda;
558 m_B[0][0][0] = 1.0 - lambda;
559 m_B[0][0][1] = lambda;
573 m_A = Array<OneD, Array<TwoD,NekDouble> >(1);
574 m_B = Array<OneD, Array<TwoD,NekDouble> >(1);
583 m_A[0][0][0] = lambda;
584 m_A[0][1][0] = 0.5 * (1.0 - lambda);
585 m_A[0][2][0] = 0.25 * (-6.0*lambda*lambda + 16.0*lambda - 1.0);
586 m_A[0][1][1] = lambda;
587 m_A[0][2][1] = 0.25 * ( 6.0*lambda*lambda - 20.0*lambda + 5.0);
588 m_A[0][2][2] = lambda;
590 m_B[0][0][0] = 0.25 * (-6.0*lambda*lambda + 16.0*lambda - 1.0);
591 m_B[0][0][1] = 0.25 * ( 6.0*lambda*lambda - 20.0*lambda + 5.0);
592 m_B[0][0][2] = lambda;
606 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
607 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
619 m_A[0][1][1] = lambda;
620 m_A[0][2][1] = 1.0 - lambda;
621 m_A[0][2][2] = lambda;
623 m_B[0][0][1] = 1.0 - lambda;
624 m_B[0][0][2] = lambda;
626 m_A[1][1][0] = lambda;
627 m_A[1][2][0] = delta;
628 m_A[1][2][1] = 1.0 - delta;
630 m_B[1][0][1] = 1.0 - lambda;
631 m_B[1][0][2] = lambda;
645 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
646 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
657 m_A[0][1][1] = lambda;
658 m_A[0][2][1] = 0.5 * (1.0 - lambda);
659 m_A[0][3][1] = 0.25 * (-6.0*lambda*lambda + 16.0*lambda - 1.0);
660 m_A[0][2][2] = lambda;
661 m_A[0][3][2] = 0.25 * ( 6.0*lambda*lambda - 20.0*lambda + 5.0);
662 m_A[0][3][3] = lambda;
664 m_B[0][0][1] = 0.25 * (-6.0*lambda*lambda + 16.0*lambda - 1.0);
665 m_B[0][0][2] = 0.25 * ( 6.0*lambda*lambda - 20.0*lambda + 5.0);
666 m_B[0][0][3] = lambda;
668 m_A[1][1][0] = 0.4358665215;
669 m_A[1][2][0] = 0.3212788860;
670 m_A[1][2][1] = 0.3966543747;
671 m_A[1][3][0] =-0.105858296;
672 m_A[1][3][1] = 0.5529291479;
673 m_A[1][3][2] = 0.5529291479;
675 m_B[1][0][1] = 0.25 * (-6.0*lambda*lambda + 16.0*lambda - 1.0);
676 m_B[1][0][2] = 0.25 * ( 6.0*lambda*lambda - 20.0*lambda + 5.0);
677 m_B[1][0][3] = lambda;
692 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
693 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
702 m_B[0][0][0] = secondth;
705 m_U[0][0] = 2*secondth;
706 m_U[0][2] = 3*secondth;
707 m_U[0][3] = -1*secondth;
709 m_V[0][0] = 2*secondth;
710 m_V[0][1] = secondth;
711 m_V[0][2] = 3*secondth;
712 m_V[0][3] = -1*secondth;
720 m_timeLevelOffset[1] = 0;
721 m_timeLevelOffset[2] = 0;
722 m_timeLevelOffset[3] = 1;
731 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
732 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
741 m_B[0][0][0] = twothirdth;
742 m_B[1][0][0] = twothirdth;
743 m_U[0][0] = 2*twothirdth;
744 m_U[0][1] = -0.5*twothirdth;
746 m_V[0][0] = 2*twothirdth;
747 m_V[0][1] = -0.5*twothirdth;
755 m_timeLevelOffset[1] = 1;
764 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
765 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
774 m_B[0][0][0] = sixthx;
778 m_U[0][1] = 6.0/16.0;
779 m_U[0][2] = 1.0/16.0;
784 m_V[0][1] = 6.0/16.0;
785 m_V[0][2] = 1.0/16.0;
796 m_timeLevelOffset[1] = 0;
797 m_timeLevelOffset[2] = 1;
798 m_timeLevelOffset[3] = 0;
799 m_timeLevelOffset[4] = 1;
807 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
808 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
820 m_A[0][1][1] = glambda;
821 m_A[0][2][1] = 1.0 - glambda;
822 m_A[0][2][2] = glambda;
824 m_B[0][0][1] = 1.0 - glambda;
825 m_B[0][0][2] = glambda;
827 m_A[1][1][0] = glambda;
828 m_A[1][2][0] = gdelta;
829 m_A[1][2][1] = 1.0 - gdelta;
831 m_B[1][0][0] = gdelta;
832 m_B[1][0][1] = 1.0 - gdelta;
846 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
847 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
858 m_A[0][1][1] = glambda;
859 m_A[0][2][1] = 1.0 - 2.0*glambda;
860 m_A[0][2][2] = glambda;
865 m_A[1][1][0] = glambda;
866 m_A[1][2][0] = glambda - 1.0;
867 m_A[1][2][1] = 2.0*(1-glambda);
884 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
885 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
914 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
915 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
945 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
946 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
955 m_A[0][1][1] = 1.0/2.0;
959 m_A[1][1][0] = 1.0/2.0;
975 m_A = Array<OneD, Array<TwoD,NekDouble> >(2);
976 m_B = Array<OneD, Array<TwoD,NekDouble> >(2);
985 m_A[0][1][1] = 1.0/2.0;
986 m_A[0][2][1] = 1.0/6.0;
987 m_A[0][2][2] = 1.0/2.0;
988 m_A[0][3][1] = -1.0/2.0;
989 m_A[0][3][2] = 1.0/2.0;
990 m_A[0][3][3] = 1.0/2.0;
991 m_A[0][4][1] = 3.0/2.0;
992 m_A[0][4][2] = -3.0/2.0;
993 m_A[0][4][3] = 1.0/2.0;
994 m_A[0][4][4] = 1.0/2.0;
997 m_B[0][0][1] = 3.0/2.0;
998 m_B[0][0][2] = -3.0/2.0;
999 m_B[0][0][3] = 1.0/2.0;
1000 m_B[0][0][4] = 1.0/2.0;
1002 m_A[1][1][0] = 1.0/2.0;
1003 m_A[1][2][0] = 11.0/18.0;
1004 m_A[1][2][1] = 1.0/18.0;
1005 m_A[1][3][0] = 5.0/6.0;
1006 m_A[1][3][1] = -5.0/6.0;
1007 m_A[1][3][2] = 1.0/2.0;
1008 m_A[1][4][0] = 1.0/4.0;
1009 m_A[1][4][1] = 7.0/4.0;
1010 m_A[1][4][2] = 3.0/4.0;
1011 m_A[1][4][3] = -7.0/4.0;
1013 m_B[1][0][0] = 1.0/4.0;
1014 m_B[1][0][1] = 7.0/4.0;
1015 m_B[1][0][2] = 3.0/4.0;
1016 m_B[1][0][3] = -7.0/4.0;
1035 "Time integration scheme coefficients do not match its type");
1041 const Array<
OneD,
const Array<TwoD, NekDouble> >& A,
1042 const Array<
OneD,
const Array<TwoD, NekDouble> >& B,
1043 const Array<TwoD, const NekDouble>& U,
1044 const Array<TwoD, const NekDouble>& V)
const
1049 int IMEXdim =
A.num_elements();
1050 int dim =
A[0].GetRows();
1052 Array<OneD, TimeIntegrationSchemeType> vertype(IMEXdim,
eExplicit);
1054 for(m = 0; m < IMEXdim; m++)
1056 for(i = 0; i < dim; i++)
1064 for(i = 0; i < dim; i++)
1066 for(j = i+1; j < dim; j++)
1071 ASSERTL1(
false,
"Fully Impplicit schemes cannnot be handled by the TimeIntegrationScheme class");
1079 ASSERTL1(
B.num_elements()==2,
"Coefficient Matrix B should have an implicit and explicit part for IMEX schemes");
1087 ASSERTL1(
false,
"This is not a proper IMEX scheme");
1091 return (vertype[0] == type);
1113 int nvar = y_0.num_elements();
1114 int npoints = y_0[0].num_elements();
1116 for(i = 0; i < nvar; i++)
1118 f_y_0[i] = Array<OneD,NekDouble>(npoints);
1124 for(i = 0; i < nvar; i++)
1126 Blas::Dscal(npoints,timestep,f_y_0[i].
get(),1);
1128 y_out->SetDerivative(0,f_y_0,timestep);
1141 "Fully Implicit integration scheme cannot be handled by this routine.");
1143 int nvar = solvector->GetFirstDim ();
1144 int npoints = solvector->GetSecondDim();
1146 if( (solvector->GetIntegrationScheme()).
get() != this )
1182 unsigned int nMasterSchemeVals = solvector->GetNvalues();
1183 unsigned int nMasterSchemeDers = solvector->GetNderivs();
1188 const Array<OneD, const unsigned int>& masterTimeLevels = solvector->GetTimeLevelOffset();
1196 for(n = 0; n < nCurSchemeVals; n++)
1202 y_n = solvector->GetValue ( curTimeLevels[n] );
1203 t_n = solvector->GetValueTime( curTimeLevels[n] );
1207 solvector_in->SetValue(curTimeLevels[n],y_n,t_n);
1209 for(n = nCurSchemeVals; n < nCurSchemeSteps; n++)
1214 dtFy_n = solvector->GetDerivative ( curTimeLevels[n] );
1218 solvector_in->SetDerivative(curTimeLevels[n],dtFy_n,timestep);
1227 solvector_in->GetTimeVector(),
1228 solvector_out->UpdateSolutionVector(),
1229 solvector_out->UpdateTimeVector(),op);
1243 bool CalcNewDeriv =
false;
1245 if( nMasterSchemeDers > 0 )
1247 if(nCurSchemeDers == 0)
1249 CalcNewDeriv =
true;
1253 if( masterTimeLevels[nMasterSchemeVals] < curTimeLevels[nCurSchemeVals] )
1255 CalcNewDeriv =
true;
1262 int newDerivTimeLevel = masterTimeLevels[nMasterSchemeVals];
1269 if (newDerivTimeLevel == 0)
1271 y_n = solvector_out->GetValue(0);
1272 t_n = solvector_out->GetValueTime(0);
1276 else if( newDerivTimeLevel == 1 )
1278 y_n = solvector->GetValue(0);
1279 t_n = solvector->GetValueTime(0);
1283 ASSERTL1(
false,
"Problems with initialising scheme");
1287 for(j = 0; j < nvar; j++)
1289 f_n[j] = Array<OneD, NekDouble>(npoints);
1296 for(j = 0; j < nvar; j++)
1304 solvector->RotateSolutionVector();
1306 solvector->SetDerivative(newDerivTimeLevel,f_n,timestep);
1312 solvector->RotateSolutionVector();
1319 for(n = 0; n < nCurSchemeVals; n++)
1325 y_n = solvector_out->GetValue ( curTimeLevels[n] );
1326 t_n = solvector_out->GetValueTime( curTimeLevels[n] );
1329 solvector->SetValue(curTimeLevels[n],y_n,t_n);
1332 for(n = nCurSchemeVals; n < nCurSchemeSteps; n++)
1338 dtFy_n = solvector_out->GetDerivative ( curTimeLevels[n] );
1342 solvector->SetDerivative(curTimeLevels[n],dtFy_n,timestep);
1352 solvector->GetTimeVector(),
1353 solvector_new->UpdateSolutionVector(),
1354 solvector_new->UpdateTimeVector(),op);
1356 solvector = solvector_new;
1358 return solvector->GetSolution();
1378 for(j = 0; j <
m_nvar; j++)
1402 for(j = 0; j <
m_nvar; j++)
1416 for(j = 0; j <
m_nvar; j++)
1428 for(j = 0; j <
m_nvar; j++)
1440 for(j = 0; j <
m_nvar; j++)
1456 for(k = 0; k <
m_nvar; k++)
1469 for(k = 0; k <
m_nvar; k++)
1482 m_T =
A(i,0)*timestep;
1484 for( j = 1; j < i; j++ )
1486 for(k = 0; k <
m_nvar; k++)
1498 m_T +=
A(i,j)*timestep;
1505 for(k = 0; k <
m_nvar; k++)
1510 m_T +=
U(i,j)*t_old[j];
1519 m_T= t_old[0]+timestep;
1524 for(
int j=0; j<=i; ++j)
1526 m_T +=
A(i,j)*timestep;
1532 for(k = 0; k <
m_nvar; k++)
1538 else if(type ==
eIMEX)
1542 m_T= t_old[0]+timestep;
1547 for(
int j=0; j<=i; ++j)
1549 m_T +=
A(i,j)*timestep;
1557 for(k = 0; k <
m_nvar; k++)
1561 m_F[i][k],1,
m_F[i][k],1);
1589 for(k = 0; k <
m_nvar; k++)
1594 if (m_numstages==1 && type ==
eIMEX)
1596 t_new[0] = t_old[0]+timestep;
1600 t_new[0] =
B(0,0)*timestep;
1603 t_new[0] +=
B(0,j)*timestep;
1607 t_new[0] +=
V(0,j)*t_old[j];
1618 for(k = 0; k <
m_nvar; k++)
1630 if(m_numstages != 1 || type !=
eIMEX)
1632 t_new[i] =
B(i,0)*timestep;
1638 for(k = 0; k <
m_nvar; k++)
1641 y_new[i][k],1,y_new[i][k],1);
1650 if(m_numstages != 1 || type !=
eIMEX)
1652 t_new[i] +=
B(i,j)*timestep;
1660 for(k = 0; k <
m_nvar; k++)
1663 y_new[i][k],1,y_new[i][k],1);
1665 if(m_numstages != 1 || type !=
eIMEX)
1667 t_new[i] +=
V(i,j)*t_old[j];
1680 const Array<
OneD,
const Array<TwoD, NekDouble> >& B,
1681 const Array<TwoD, const NekDouble>& U,
1682 const Array<TwoD, const NekDouble>& V)
const
1690 for(m = 0; m <
A.num_elements(); m++)
1718 const Array<
OneD,
const Array<TwoD, NekDouble> >& B,
1719 const Array<TwoD, const NekDouble>& U,
1720 const Array<TwoD, const NekDouble>& V)
const
1728 for(m = 0; m <
A.num_elements(); m++)
1762 ASSERTL1(y_old[0]. num_elements()==y_new[0]. num_elements(),
"Non-matching number of variables.");
1763 ASSERTL1(y_old[0][0].num_elements()==y_new[0][0].num_elements(),
"Non-matching number of coefficients.");
1784 int osprecision = 6;
1787 os <<
"- number of steps: " << r << endl;
1788 os <<
"- number of stages: " << s << endl;
1790 os <<
"General linear method tableau: " << endl;
1792 for(i = 0; i < s; i++)
1794 for(j = 0; j < s; j++)
1797 os.precision(osprecision);
1798 os << right << rhs.
A(i,j) <<
" ";
1803 for(j = 0; j < s; j++)
1806 os.precision(osprecision);
1807 os << right << rhs.
A_IMEX(i,j) <<
" ";
1812 for(j = 0; j < r; j++)
1815 os.precision(osprecision);
1816 os << right << rhs.
U(i,j);
1820 int imexflag = (type ==
eIMEX)?2:1;
1821 for(
int i = 0; i < (r+imexflag*s)*(oswidth+1)+imexflag*2-1; i++)
1826 for(i = 0; i < r; i++)
1828 for(j = 0; j < s; j++)
1831 os.precision(osprecision);
1832 os << right << rhs.
B(i,j) <<
" ";
1837 for(j = 0; j < s; j++)
1840 os.precision(osprecision);
1841 os << right << rhs.
B_IMEX(i,j) <<
" ";
1846 for(j = 0; j < r; j++)
1849 os.precision(osprecision);
1850 os << right << rhs.
V(i,j);