void GatherCoordinateVector<PHAL::AlbanyTraits::Tangent, Traits>:: evaluateFields(typename Traits::EvalData workset) { unsigned int numCells = workset.numCells; Teuchos::ArrayRCP<Teuchos::ArrayRCP<double*> > wsCoords = workset.wsCoords; const Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::ArrayRCP<double> > > > & ws_coord_derivs = workset.ws_coord_derivs; std::vector<int> *coord_deriv_indices = workset.coord_deriv_indices; int numShapeDerivs = ws_coord_derivs.size(); int numParams = workset.num_cols_p; for (std::size_t cell=0; cell < numCells; ++cell) { for (std::size_t node = 0; node < numVertices; ++node) { for (std::size_t eq=0; eq < numDim; ++eq) { coordVec(cell,node,eq) = TanFadType(numParams, wsCoords[cell][node][eq]); for (int j=0; j < numShapeDerivs; ++j) { coordVec(cell,node,eq).fastAccessDx((*coord_deriv_indices)[j]) = ws_coord_derivs[j][cell][node][eq]; } } } } // Since Intrepid2 will later perform calculations on the entire workset size // and not just the used portion, we must fill the excess with reasonable // values. Leaving this out leads to calculations on singular elements. // for (std::size_t cell=numCells; cell < worksetSize; ++cell) { for (std::size_t node = 0; node < numVertices; ++node) { for (std::size_t eq=0; eq < numDim; ++eq) { coordVec(cell,node,eq) = coordVec(0,node,eq); } } } }
void GatherSolution<PHAL::AlbanyTraits::Tangent, Traits>:: evaluateFields(typename Traits::EvalData workset) { Teuchos::RCP<const Epetra_Vector> x = workset.x; Teuchos::RCP<const Epetra_Vector> xdot = workset.xdot; Teuchos::RCP<const Epetra_MultiVector> Vx = workset.Vx; Teuchos::RCP<const Epetra_MultiVector> Vxdot = workset.Vxdot; Teuchos::RCP<ParamVec> params = workset.params; int num_cols_tot = workset.param_offset + workset.num_cols_p; ScalarT* valptr; for (int cell=0; cell < workset.numCells; ++cell ) { const Teuchos::ArrayRCP<Teuchos::ArrayRCP<int> >& nodeID = workset.wsElNodeEqID[cell]; for (int node = 0; node < this->numNodes; ++node) { const Teuchos::ArrayRCP<int>& eqID = nodeID[node]; int n = 0, eq = 0; for (int j = eq; j < eq+this->numNodeVar; j++, ++n) { valptr = &(this->val[j])(cell,node); if (Vx != Teuchos::null && workset.j_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*x)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.j_coeff*(*Vx)[k][eqID[n]]; } else *valptr = TanFadType((*x)[eqID[n]]); } eq += this->numNodeVar; for (int level = 0; level < this->numLevels; level++) { for (int j = eq; j < eq+this->numVectorLevelVar; j++) { for (int dim = 0; dim < this->numDims; ++dim, ++n) { valptr = &(this->val[j])(cell,node,level,dim); if (Vx != Teuchos::null && workset.j_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*x)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.j_coeff*(*Vx)[k][eqID[n]]; } else *valptr = TanFadType((*x)[eqID[n]]); } } for (int j = eq+this->numVectorLevelVar; j < eq+this->numVectorLevelVar+this->numScalarLevelVar; j++, ++n) { valptr = &(this->val[j])(cell,node,level); if (Vx != Teuchos::null && workset.j_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*x)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.j_coeff*(*Vx)[k][eqID[n]]; } else *valptr = TanFadType((*x)[eqID[n]]); } } eq += this->numVectorLevelVar+this->numScalarLevelVar; for (int level = 0; level < this->numLevels; ++level) { for (int j = eq; j < eq+this->numTracerVar; ++j, ++n) { valptr = &(this->val[j])(cell,node,level); if (Vx != Teuchos::null && workset.j_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*x)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.j_coeff*(*Vx)[k][eqID[n]]; } else *valptr = TanFadType((*x)[eqID[n]]); } } eq += this->numTracerVar; if (workset.transientTerms) { int n = 0, eq = 0; for (int j = eq; j < eq+this->numNodeVar; j++, ++n) { valptr = &(this->val_dot[j])(cell,node); if (Vxdot != Teuchos::null && workset.m_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*xdot)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.m_coeff*(*Vxdot)[k][eqID[n]]; } else *valptr = TanFadType((*xdot)[eqID[n]]); } eq += this->numNodeVar; for (int level = 0; level < this->numLevels; level++) { for (int j = eq; j < eq+this->numVectorLevelVar; j++) { for (int dim = 0; dim < this->numDims; ++dim, ++n) { valptr = &(this->val_dot[j])(cell,node,level,dim); if (Vxdot != Teuchos::null && workset.m_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*xdot)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.m_coeff*(*Vxdot)[k][eqID[n]]; } else *valptr = TanFadType((*xdot)[eqID[n]]); } } for (int j = eq+this->numVectorLevelVar; j < eq+this->numScalarLevelVar+this->numScalarLevelVar; j++,++n) { valptr = &(this->val_dot[j])(cell,node,level); if (Vxdot != Teuchos::null && workset.m_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*xdot)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.m_coeff*(*Vxdot)[k][eqID[n]]; } else *valptr = TanFadType((*xdot)[eqID[n]]); } } eq += this->numVectorLevelVar+this->numScalarLevelVar; for (int level = 0; level < this->numLevels; ++level) { for (int j = eq; j < eq+this->numTracerVar; ++j, ++n) { valptr = &(this->val_dot[j])(cell,node,level); if (Vxdot != Teuchos::null && workset.m_coeff != 0.0) { *valptr = TanFadType(num_cols_tot, (*xdot)[eqID[n]]); for (int k=0; k<workset.num_cols_x; k++) valptr->fastAccessDx(k) = workset.m_coeff*(*Vxdot)[k][eqID[n]]; } else *valptr = TanFadType((*xdot)[eqID[n]]); } } eq += this->numTracerVar; } } } }