72 m_substeps = pSession->GetParameter(
"Substeps");
78 std::set<enum LibUtilities::ShapeType> s;
79 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
81 s.insert(m_field->GetExp(i)->DetShapeType());
96 m_nq = pField->GetNcoeffs();
97 int order = m_field->GetExp(0)->GetBasis(0)->GetNumModes();
117 m_nq = pField->GetTotPoints();
125 void CellModel::Initialise()
127 ASSERTL1(m_nvar > 0,
"Cell model must have at least 1 variable.");
131 for (
unsigned int i = 0; i < m_nvar; ++i)
137 for (
unsigned int i = 0; i < m_gates.size(); ++i)
142 if (m_session->DefinesFunction(
"CellModelInitialConditions"))
148 v_SetInitialConditions();
160 void CellModel::TimeIntegrate(
167 int nvar = inarray.size();
174 if (!m_nodalTmp.size())
177 for (
unsigned int k = 0; k < nvar; ++k)
186 for (
unsigned int k = 0; k < nvar; ++k)
188 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
190 phys_offset = m_field->GetPhys_Offset(i);
191 coef_offset = m_field->GetCoeff_Offset(i);
194 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
195 m_nodalTri->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
199 m_field->GetExp(0)->FwdTrans(inarray[k] + phys_offset, tmpCoeffs);
200 m_nodalTet->ModalToNodal(tmpCoeffs, tmp=m_nodalTmp[k]+coef_offset);
214 NekDouble delta_t = (time - m_lastTime)/m_substeps;
218 for (
unsigned int i = 0; i < m_substeps - 1; ++i)
220 Update(m_cellSol, m_wsp, time);
222 Vmath::Svtvp(m_nq, delta_t, m_wsp[0], 1, m_cellSol[0], 1, m_cellSol[0], 1);
224 for (
unsigned int j = 0; j < m_concentrations.size(); ++j)
226 Vmath::Svtvp(m_nq, delta_t, m_wsp[m_concentrations[j]], 1, m_cellSol[m_concentrations[j]], 1, m_cellSol[m_concentrations[j]], 1);
229 for (
unsigned int j = 0; j < m_gates.size(); ++j)
231 Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
232 Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
233 Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
234 Vmath::Vvtvp(m_nq, m_cellSol[m_gates[j]], 1, m_gates_tau[j], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
239 Update(m_cellSol, m_wsp, time);
247 for (
unsigned int k = 0; k < nvar; ++k)
249 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
251 int phys_offset = m_field->GetPhys_Offset(i);
252 int coef_offset = m_field->GetCoeff_Offset(i);
255 m_nodalTri->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
256 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
260 m_nodalTet->NodalToModal(m_wsp[k]+coef_offset, tmpCoeffs);
261 m_field->GetExp(0)->BwdTrans(tmpCoeffs, tmp=outarray[k] + phys_offset);
272 for (
unsigned int j = 0; j < m_concentrations.size(); ++j)
274 Vmath::Svtvp(m_nq, delta_t, m_wsp[m_concentrations[j]], 1, m_cellSol[m_concentrations[j]], 1, m_cellSol[m_concentrations[j]], 1);
278 for (
unsigned int j = 0; j < m_gates.size(); ++j)
280 Vmath::Sdiv(m_nq, -delta_t, m_gates_tau[j], 1, m_gates_tau[j], 1);
281 Vmath::Vexp(m_nq, m_gates_tau[j], 1, m_gates_tau[j], 1);
282 Vmath::Vsub(m_nq, m_cellSol[m_gates[j]], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
283 Vmath::Vvtvp(m_nq, m_cellSol[m_gates[j]], 1, m_gates_tau[j], 1, m_wsp[m_gates[j]], 1, m_cellSol[m_gates[j]], 1);
291 ASSERTL0(idx < m_nvar,
"Index out of range for cell model.");
298 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
300 int coef_offset = m_field->GetCoeff_Offset(i);
303 m_nodalTri->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
307 m_nodalTet->NodalToModal(m_cellSol[idx]+coef_offset, tmp=outarray+coef_offset);
313 m_field->FwdTrans_IterPerExp(m_cellSol[idx], outarray);
321 return m_cellSol[idx];
324 void CellModel::LoadCellModel()
326 const bool root = (m_session->GetComm()->GetRank() == 0);
327 const std::string fncName =
"CellModelInitialConditions";
328 const int nvar = m_cellSol[0].size();
338 cout <<
"Cell model initial conditions: " << endl;
342 std::set<std::string> filelist;
343 for (j = 1; j < nvar; ++j)
346 varName = GetCellVarName(j);
348 if (m_session->GetFunctionType(fncName, varName) ==
351 filelist.insert(m_session->GetFunctionFilename(fncName,
357 typedef std::vector<LibUtilities::FieldDefinitionsSharedPtr> FDef;
358 typedef std::vector<std::vector<NekDouble> > FData;
359 std::map<std::string, FDef> FieldDef;
360 std::map<std::string, FData> FieldData;
363 for (
auto &setIt : filelist)
367 cout <<
" - Reading file: " << setIt << endl;
369 FieldDef[setIt] = FDef(0);
370 FieldData[setIt] = FData(0);
372 LibUtilities::FieldIO::CreateForFile(m_session, setIt);
373 fld->Import(setIt, FieldDef[setIt], FieldData[setIt],
378 auto iter = fieldMetaDataMap.find(
"Time");
379 if(iter != fieldMetaDataMap.end())
381 m_lastTime = boost::lexical_cast<NekDouble>(iter->second);
387 for(j = 1; j < m_cellSol.size(); ++j)
390 varName = GetCellVarName(j);
393 if (m_session->GetFunctionType(fncName, varName) ==
396 const std::string file =
397 m_session->GetFunctionFilename(fncName, varName);
401 cout <<
" - Field " << varName <<
": from file "
406 for(
int i = 0; i < FieldDef[file].size(); ++i)
408 m_field->ExtractDataToCoeffs(FieldDef[file][i],
418 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
420 int coef_offset = m_field->GetCoeff_Offset(i);
421 if (m_field->GetExp(0)->DetShapeType() ==
424 m_nodalTri->ModalToNodal(coeffs+coef_offset,
425 tmp=m_cellSol[j]+coef_offset);
429 m_nodalTet->ModalToNodal(coeffs+coef_offset,
430 tmp=m_cellSol[j]+coef_offset);
436 m_field->BwdTrans(coeffs, m_cellSol[j]);
439 else if (m_session->GetFunctionType(fncName, varName) ==
443 m_session->GetFunction(fncName, varName);
447 cout <<
" - Field " << varName <<
": "
448 << equ->GetExpression() << endl;
451 const unsigned int nphys = m_field->GetNpoints();
455 m_field->GetCoords(x0,x1,x2);
462 equ->Evaluate(x0, x1, x2, phys);
463 for (
unsigned int i = 0; i < m_field->GetNumElmts(); ++i)
465 int phys_offset = m_field->GetPhys_Offset(i);
466 int coef_offset = m_field->GetCoeff_Offset(i);
467 if (m_field->GetExp(0)->DetShapeType() ==
470 m_field->GetExp(0)->FwdTrans(
471 phys + phys_offset, tmpCoeffs);
472 m_nodalTri->ModalToNodal(
474 tmp = m_cellSol[j] + coef_offset);
478 m_field->GetExp(0)->FwdTrans(
479 phys + phys_offset, tmpCoeffs);
480 m_nodalTet->ModalToNodal(
482 tmp = m_cellSol[j] + coef_offset);
488 equ->Evaluate(x0, x1, x2, m_cellSol[j]);
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
Provides a generic Factory class.
Defines a specification for a set of points.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
std::shared_ptr< FieldIO > FieldIOSharedPtr
std::map< std::string, std::string > FieldMetaDataMap
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
@ eNodalTriEvenlySpaced
2D Evenly-spaced points on a Triangle
@ eGaussRadauMAlpha2Beta0
Gauss Radau pinned at x=-1, .
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ eNodalTetEvenlySpaced
3D Evenly-spaced points on a Tetrahedron
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
@ eFunctionTypeExpression
@ eModified_B
Principle Modified Functions .
@ eModified_C
Principle Modified Functions .
@ eModified_A
Principle Modified Functions .
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
CellModelFactory & GetCellModelFactory()
void Vexp(int n, const T *x, const int incx, T *y, const int incy)
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.