void AssembleBoussinesqAppoximation_AD(MultiLevelProblem& ml_prob) { // ml_prob is the global object from/to where get/set all the data // level is the level of the PDE system to be assembled // levelMax is the Maximum level of the MultiLevelProblem // assembleMatrix is a flag that tells if only the residual or also the matrix should be assembled // extract pointers to the several objects that we are going to use TransientNonlinearImplicitSystem* mlPdeSys = &ml_prob.get_system<TransientNonlinearImplicitSystem> ("NS"); // pointer to the linear implicit system named "Poisson" const unsigned level = mlPdeSys->GetLevelToAssemble(); Mesh* msh = ml_prob._ml_msh->GetLevel(level); // pointer to the mesh (level) object elem* el = msh->el; // pointer to the elem object in msh (level) MultiLevelSolution* mlSol = ml_prob._ml_sol; // pointer to the multilevel solution object Solution* sol = ml_prob._ml_sol->GetSolutionLevel(level); // pointer to the solution (level) object LinearEquationSolver* pdeSys = mlPdeSys->_LinSolver[level]; // pointer to the equation (level) object bool assembleMatrix = mlPdeSys->GetAssembleMatrix(); // call the adept stack object adept::Stack& s = FemusInit::_adeptStack; if(assembleMatrix) s.continue_recording(); else s.pause_recording(); SparseMatrix* KK = pdeSys->_KK; // pointer to the global stifness matrix object in pdeSys (level) NumericVector* RES = pdeSys->_RES; // pointer to the global residual vector object in pdeSys (level) const unsigned dim = msh->GetDimension(); // get the domain dimension of the problem unsigned dim2 = (3 * (dim - 1) + !(dim - 1)); // dim2 is the number of second order partial derivatives (1,3,6 depending on the dimension) unsigned iproc = msh->processor_id(); // get the process_id (for parallel computation) // reserve memory for the local standar vectors const unsigned maxSize = static_cast< unsigned >(ceil(pow(3, dim))); // conservative: based on line3, quad9, hex27 //solution variable unsigned solTIndex; solTIndex = mlSol->GetIndex("T"); // get the position of "T" in the ml_sol object unsigned solTType = mlSol->GetSolutionType(solTIndex); // get the finite element type for "T" vector < unsigned > solVIndex(dim); solVIndex[0] = mlSol->GetIndex("U"); // get the position of "U" in the ml_sol object solVIndex[1] = mlSol->GetIndex("V"); // get the position of "V" in the ml_sol object if(dim == 3) solVIndex[2] = mlSol->GetIndex("W"); // get the position of "V" in the ml_sol object unsigned solVType = mlSol->GetSolutionType(solVIndex[0]); // get the finite element type for "u" unsigned solPIndex; solPIndex = mlSol->GetIndex("P"); // get the position of "P" in the ml_sol object unsigned solPType = mlSol->GetSolutionType(solPIndex); // get the finite element type for "u" unsigned solTPdeIndex; solTPdeIndex = mlPdeSys->GetSolPdeIndex("T"); // get the position of "T" in the pdeSys object // std::cout << solTIndex <<" "<<solTPdeIndex<<std::endl; vector < unsigned > solVPdeIndex(dim); solVPdeIndex[0] = mlPdeSys->GetSolPdeIndex("U"); // get the position of "U" in the pdeSys object solVPdeIndex[1] = mlPdeSys->GetSolPdeIndex("V"); // get the position of "V" in the pdeSys object if(dim == 3) solVPdeIndex[2] = mlPdeSys->GetSolPdeIndex("W"); unsigned solPPdeIndex; solPPdeIndex = mlPdeSys->GetSolPdeIndex("P"); // get the position of "P" in the pdeSys object vector < adept::adouble > solT; // local solution vector < vector < adept::adouble > > solV(dim); // local solution vector < adept::adouble > solP; // local solution vector < double > solTold; // local solution vector < vector < double > > solVold(dim); // local solution vector < double > solPold; // local solution vector< adept::adouble > aResT; // local redidual vector vector< vector < adept::adouble > > aResV(dim); // local redidual vector vector< adept::adouble > aResP; // local redidual vector vector < vector < double > > coordX(dim); // local coordinates unsigned coordXType = 2; // get the finite element type for "x", it is always 2 (LAGRANGE QUADRATIC) solT.reserve(maxSize); solTold.reserve(maxSize); aResT.reserve(maxSize); for(unsigned k = 0; k < dim; k++) { solV[k].reserve(maxSize); solVold[k].reserve(maxSize); aResV[k].reserve(maxSize); coordX[k].reserve(maxSize); } solP.reserve(maxSize); solPold.reserve(maxSize); aResP.reserve(maxSize); vector <double> phiV; // local test function vector <double> phiV_x; // local test function first order partial derivatives vector <double> phiV_xx; // local test function second order partial derivatives phiV.reserve(maxSize); phiV_x.reserve(maxSize * dim); phiV_xx.reserve(maxSize * dim2); vector <double> phiT; // local test function vector <double> phiT_x; // local test function first order partial derivatives vector <double> phiT_xx; // local test function second order partial derivatives phiT.reserve(maxSize); phiT_x.reserve(maxSize * dim); phiT_xx.reserve(maxSize * dim2); double* phiP; double weight; // gauss point weight vector< int > sysDof; // local to global pdeSys dofs sysDof.reserve((dim + 2) *maxSize); vector< double > Res; // local redidual vector Res.reserve((dim + 2) *maxSize); vector < double > Jac; Jac.reserve((dim + 2) *maxSize * (dim + 2) *maxSize); if(assembleMatrix) KK->zero(); // Set to zero all the entries of the Global Matrix // element loop: each process loops only on the elements that owns for(int iel = msh->_elementOffset[iproc]; iel < msh->_elementOffset[iproc + 1]; iel++) { // element geometry type short unsigned ielGeom = msh->GetElementType(iel); unsigned nDofsT = msh->GetElementDofNumber(iel, solTType); // number of solution element dofs unsigned nDofsV = msh->GetElementDofNumber(iel, solVType); // number of solution element dofs unsigned nDofsP = msh->GetElementDofNumber(iel, solPType); // number of solution element dofs unsigned nDofsX = msh->GetElementDofNumber(iel, coordXType); // number of coordinate element dofs unsigned nDofsTVP = nDofsT + dim * nDofsV + nDofsP; // resize local arrays sysDof.resize(nDofsTVP); solT.resize(nDofsV); solTold.resize(nDofsV); for(unsigned k = 0; k < dim; k++) { solV[k].resize(nDofsV); solVold[k].resize(nDofsV); coordX[k].resize(nDofsX); } solP.resize(nDofsP); solPold.resize(nDofsP); aResT.resize(nDofsV); //resize std::fill(aResT.begin(), aResT.end(), 0); //set aRes to zero for(unsigned k = 0; k < dim; k++) { aResV[k].resize(nDofsV); //resize std::fill(aResV[k].begin(), aResV[k].end(), 0); //set aRes to zero } aResP.resize(nDofsP); //resize std::fill(aResP.begin(), aResP.end(), 0); //set aRes to zero // local storage of global mapping and solution for(unsigned i = 0; i < nDofsT; i++) { unsigned solTDof = msh->GetSolutionDof(i, iel, solTType); // global to global mapping between solution node and solution dof solT[i] = (*sol->_Sol[solTIndex])(solTDof); // global extraction and local storage for the solution solTold[i] = (*sol->_SolOld[solTIndex])(solTDof); // global extraction and local storage for the solution sysDof[i] = pdeSys->GetSystemDof(solTIndex, solTPdeIndex, i, iel); // global to global mapping between solution node and pdeSys dofs } // local storage of global mapping and solution for(unsigned i = 0; i < nDofsV; i++) { unsigned solVDof = msh->GetSolutionDof(i, iel, solVType); // global to global mapping between solution node and solution dof for(unsigned k = 0; k < dim; k++) { solV[k][i] = (*sol->_Sol[solVIndex[k]])(solVDof); // global extraction and local storage for the solution solVold[k][i] = (*sol->_SolOld[solVIndex[k]])(solVDof); // global extraction and local storage for the solution sysDof[i + nDofsT + k * nDofsV] = pdeSys->GetSystemDof(solVIndex[k], solVPdeIndex[k], i, iel); // global to global mapping between solution node and pdeSys dof } } for(unsigned i = 0; i < nDofsP; i++) { unsigned solPDof = msh->GetSolutionDof(i, iel, solPType); // global to global mapping between solution node and solution dof solP[i] = (*sol->_Sol[solPIndex])(solPDof); // global extraction and local storage for the solution solPold[i] = (*sol->_SolOld[solPIndex])(solPDof); // global extraction and local storage for the solution sysDof[i + nDofsT + dim * nDofsV] = pdeSys->GetSystemDof(solPIndex, solPPdeIndex, i, iel); // global to global mapping between solution node and pdeSys dof } // local storage of coordinates for(unsigned i = 0; i < nDofsX; i++) { unsigned coordXDof = msh->GetSolutionDof(i, iel, coordXType); // global to global mapping between coordinates node and coordinate dof for(unsigned k = 0; k < dim; k++) { coordX[k][i] = (*msh->_topology->_Sol[k])(coordXDof); // global extraction and local storage for the element coordinates } } // start a new recording of all the operations involving adept::adouble variables if(assembleMatrix) s.new_recording(); // *** Gauss point loop *** for(unsigned ig = 0; ig < msh->_finiteElement[ielGeom][solVType]->GetGaussPointNumber(); ig++) { // *** get gauss point weight, test function and test function partial derivatives *** msh->_finiteElement[ielGeom][solTType]->Jacobian(coordX, ig, weight, phiT, phiT_x, phiT_xx); msh->_finiteElement[ielGeom][solVType]->Jacobian(coordX, ig, weight, phiV, phiV_x, phiV_xx); phiP = msh->_finiteElement[ielGeom][solPType]->GetPhi(ig); // evaluate the solution, the solution derivatives and the coordinates in the gauss point adept::adouble solT_gss = 0; double solTold_gss = 0; vector < adept::adouble > gradSolT_gss(dim, 0.); vector < double > gradSolTold_gss(dim, 0.); for(unsigned i = 0; i < nDofsT; i++) { solT_gss += phiT[i] * solT[i]; solTold_gss += phiT[i] * solTold[i]; for(unsigned j = 0; j < dim; j++) { gradSolT_gss[j] += phiT_x[i * dim + j] * solT[i]; gradSolTold_gss[j] += phiT_x[i * dim + j] * solTold[i]; } } vector < adept::adouble > solV_gss(dim, 0); vector < double > solVold_gss(dim, 0); vector < vector < adept::adouble > > gradSolV_gss(dim); vector < vector < double > > gradSolVold_gss(dim); for(unsigned k = 0; k < dim; k++) { gradSolV_gss[k].resize(dim); gradSolVold_gss[k].resize(dim); std::fill(gradSolV_gss[k].begin(), gradSolV_gss[k].end(), 0); std::fill(gradSolVold_gss[k].begin(), gradSolVold_gss[k].end(), 0); } for(unsigned i = 0; i < nDofsV; i++) { for(unsigned k = 0; k < dim; k++) { solV_gss[k] += phiV[i] * solV[k][i]; solVold_gss[k] += phiV[i] * solVold[k][i]; } for(unsigned j = 0; j < dim; j++) { for(unsigned k = 0; k < dim; k++) { gradSolV_gss[k][j] += phiV_x[i * dim + j] * solV[k][i]; gradSolVold_gss[k][j] += phiV_x[i * dim + j] * solVold[k][i]; } } } adept::adouble solP_gss = 0; double solPold_gss = 0; for(unsigned i = 0; i < nDofsP; i++) { solP_gss += phiP[i] * solP[i]; solPold_gss += phiP[i] * solPold[i]; } double alpha = 1.; double beta = 1.;//40000.; double Pr = 0.71; double Ra = 340000; double dt = mlPdeSys -> GetIntervalTime(); // *** phiT_i loop *** for(unsigned i = 0; i < nDofsT; i++) { adept::adouble Temp = 0.; adept::adouble TempOld = 0.; for(unsigned j = 0; j < dim; j++) { Temp += 1. / sqrt(Ra * Pr) * alpha * phiT_x[i * dim + j] * gradSolT_gss[j]; Temp += phiT[i] * (solV_gss[j] * gradSolT_gss[j]); TempOld += 1. / sqrt(Ra * Pr) * alpha * phiT_x[i * dim + j] * gradSolTold_gss[j]; TempOld += phiT[i] * (solVold_gss[j] * gradSolTold_gss[j]); } aResT[i] += (- (solT_gss - solTold_gss) * phiT[i] / dt - 0.5 * (Temp + TempOld)) * weight; } // end phiT_i loop // *** phiV_i loop *** for(unsigned i = 0; i < nDofsV; i++) { vector < adept::adouble > NSV(dim, 0.); vector < double > NSVold(dim, 0.); for(unsigned j = 0; j < dim; j++) { for(unsigned k = 0; k < dim; k++) { NSV[k] += sqrt(Pr / Ra) * phiV_x[i * dim + j] * (gradSolV_gss[k][j] + gradSolV_gss[j][k]); NSV[k] += phiV[i] * (solV_gss[j] * gradSolV_gss[k][j]); NSVold[k] += sqrt(Pr / Ra) * phiV_x[i * dim + j] * (gradSolVold_gss[k][j] + gradSolVold_gss[j][k]); NSVold[k] += phiV[i] * (solVold_gss[j] * gradSolVold_gss[k][j]); } } for(unsigned k = 0; k < dim; k++) { NSV[k] += -solP_gss * phiV_x[i * dim + k]; NSVold[k] += -solPold_gss * phiV_x[i * dim + k]; } NSV[1] += -beta * solT_gss * phiV[i]; NSVold[1] += -beta * solTold_gss * phiV[i]; for(unsigned k = 0; k < dim; k++) { aResV[k][i] += (- (solV_gss[k] - solVold_gss[k]) * phiV[i] / dt - 0.5 * (NSV[k] + NSVold[k])) * weight; } } // end phiV_i loop // *** phiP_i loop *** for(unsigned i = 0; i < nDofsP; i++) { for(int k = 0; k < dim; k++) { aResP[i] += - (gradSolV_gss[k][k]) * phiP[i] * weight; } } // end phiP_i loop } // end gauss point loop //-------------------------------------------------------------------------------------------------------- // Add the local Matrix/Vector into the global Matrix/Vector //copy the value of the adept::adoube aRes in double Res and store them in RES Res.resize(nDofsTVP); //resize for(int i = 0; i < nDofsT; i++) { Res[i] = -aResT[i].value(); } for(int i = 0; i < nDofsV; i++) { for(unsigned k = 0; k < dim; k++) { Res[ i + nDofsT + k * nDofsV ] = -aResV[k][i].value(); } } for(int i = 0; i < nDofsP; i++) { Res[ i + nDofsT + dim * nDofsV ] = -aResP[i].value(); } RES->add_vector_blocked(Res, sysDof); //Extarct and store the Jacobian if(assembleMatrix) { Jac.resize(nDofsTVP * nDofsTVP); // define the dependent variables s.dependent(&aResT[0], nDofsT); for(unsigned k = 0; k < dim; k++) { s.dependent(&aResV[k][0], nDofsV); } s.dependent(&aResP[0], nDofsP); // define the independent variables s.independent(&solT[0], nDofsT); for(unsigned k = 0; k < dim; k++) { s.independent(&solV[k][0], nDofsV); } s.independent(&solP[0], nDofsP); // get the and store jacobian matrix (row-major) s.jacobian(&Jac[0] , true); KK->add_matrix_blocked(Jac, sysDof, sysDof); s.clear_independents(); s.clear_dependents(); } } //end element loop for each process RES->close(); if(assembleMatrix) { KK->close(); } // ***************** END ASSEMBLY ******************* }
void AssembleWillmoreFlow_AD(MultiLevelProblem& ml_prob) { // ml_prob is the global object from/to where get/set all the data // level is the level of the PDE system to be assembled // levelMax is the Maximum level of the MultiLevelProblem // assembleMatrix is a flag that tells if only the residual or also the matrix should be assembled // call the adept stack object adept::Stack& s = FemusInit::_adeptStack; // extract pointers to the several objects that we are going to use TransientNonlinearImplicitSystem* mlPdeSys = &ml_prob.get_system<TransientNonlinearImplicitSystem> ("Willmore"); // pointer to the linear implicit system named "Poisson" const unsigned level = mlPdeSys->GetLevelToAssemble(); Mesh* msh = ml_prob._ml_msh->GetLevel(level); // pointer to the mesh (level) object elem* el = msh->el; // pointer to the elem object in msh (level) MultiLevelSolution* mlSol = ml_prob._ml_sol; // pointer to the multilevel solution object Solution* sol = ml_prob._ml_sol->GetSolutionLevel(level); // pointer to the solution (level) object LinearEquationSolver* pdeSys = mlPdeSys->_LinSolver[level]; // pointer to the equation (level) object SparseMatrix* KK = pdeSys->_KK; // pointer to the global stifness matrix object in pdeSys (level) NumericVector* RES = pdeSys->_RES; // pointer to the global residual vector object in pdeSys (level) const unsigned dim = msh->GetDimension(); // get the domain dimension of the problem unsigned iproc = msh->processor_id(); // get the process_id (for parallel computation) //solution variable unsigned solRIndex[3]; solRIndex[0] = mlSol->GetIndex("X"); // get the position of "X" in the ml_sol object solRIndex[1] = mlSol->GetIndex("Y"); // get the position of "Y" in the ml_sol object solRIndex[2] = mlSol->GetIndex("Z"); // get the position of "Z" in the ml_sol object unsigned solRType[3]; solRType[0] = mlSol->GetSolutionType(solRIndex[0]); // get the finite element type for "R" solRType[1] = mlSol->GetSolutionType(solRIndex[1]); // get the finite element type for "R" solRType[2] = mlSol->GetSolutionType(solRIndex[2]); // get the finite element type for "R" unsigned solRPdeIndex[3]; solRPdeIndex[0] = mlPdeSys->GetSolPdeIndex("X"); // get the position of "X" in the pdeSys object solRPdeIndex[1] = mlPdeSys->GetSolPdeIndex("Y"); // get the position of "Y" in the pdeSys object solRPdeIndex[2] = mlPdeSys->GetSolPdeIndex("Z"); // get the position of "Z" in the pdeSys object vector < adept::adouble > solR[3]; // local solution vector < double > solR_old[3]; unsigned solHIndex; solHIndex = mlSol->GetIndex("H"); // get the position of "H" in the ml_sol object unsigned solHType = mlSol->GetSolutionType(solHIndex); // get the finite element type for "H" unsigned solHPdeIndex; solHPdeIndex = mlPdeSys->GetSolPdeIndex("H"); // get the position of "H" in the pdeSys object vector < adept::adouble > solH; // local solution vector < vector < double > > x(dim); // local coordinates unsigned xType = 2; // get the finite element type for "x", it is always 2 (LAGRANGE QUADRATIC) vector< int > sysDof; // local to global pdeSys dofs vector <double> phi; // local test function vector <double> phi_x; // local test function first order partial derivatives vector <double> phi_xx; // local test function second order partial derivatives double weight; // gauss point weight vector< double > Res; // local redidual vector vector< adept::adouble > aResR[3]; // local redidual vector vector< adept::adouble > aResH; // local redidual vector // reserve memory for the local standar vectors const unsigned maxSize = static_cast< unsigned >(ceil(pow(3, dim))); // conservative: based on line3, quad9, hex27 for (int i = 0; i < 3; i++) { solR[i].reserve(maxSize); solR_old[i].reserve(maxSize); } solH.reserve(maxSize); for (unsigned i = 0; i < dim; i++) x[i].reserve(maxSize); sysDof.reserve(4 * maxSize); phi.reserve(maxSize); phi_x.reserve(maxSize * dim); unsigned dim2 = (3 * (dim - 1) + !(dim - 1)); // dim2 is the number of second order partial derivatives (1,3,6 depending on the dimension) phi_xx.reserve(maxSize * dim2); Res.reserve(4 * maxSize); for (int i = 0; i < 3; i++) { aResR[i].reserve(maxSize); } aResH.reserve(maxSize); vector < double > Jac; // local Jacobian matrix (ordered by column, adept) Jac.reserve(4 * maxSize * 4 * maxSize); KK->zero(); // Set to zero all the entries of the Global Matrix // element loop: each process loops only on the elements that owns for (int iel = msh->_elementOffset[iproc]; iel < msh->_elementOffset[iproc + 1]; iel++) { short unsigned ielGeom = msh->GetElementType(iel); unsigned nDofs = msh->GetElementDofNumber(iel, solHType); // number of solution element dofs unsigned nDofs2 = msh->GetElementDofNumber(iel, xType); // number of coordinate element dofs // resize local arrays sysDof.resize(4 * nDofs); for (int i = 0; i < 3; i++) { solR[i].resize(nDofs); solR_old[i].resize(nDofs); } solH.resize(nDofs); for (int i = 0; i < dim; i++) { x[i].resize(nDofs2); } Res.resize(4 * nDofs); //resize for (int i = 0; i < 3; i++) { aResR[i].resize(nDofs); //resize } aResH.resize(nDofs); //resize for (int i = 0; i < 3; i++) { std::fill(aResR[i].begin(), aResR[i].end(), 0); //set aRes to zero } std::fill(aResH.begin(), aResH.end(), 0); //set aRes to zero // local storage of global mapping and solution for (unsigned i = 0; i < nDofs; i++) { unsigned solDof = msh->GetSolutionDof(i, iel, solHType); // global to global mapping between solution node and solution dof for (int k = 0; k < 3; k++) { solR[k][i] = (*sol->_Sol[solRIndex[k]])(solDof); // global extraction and local storage for the solution solR_old[k][i] = (*sol->_SolOld[solRIndex[k]])(solDof); // global extraction and local storage for the solution } solH[i] = (*sol->_Sol[solHIndex])(solDof); // global extraction and local storage for the solution for (int k = 0; k < 3; k++) { sysDof[k * nDofs + i] = pdeSys->GetSystemDof(solRIndex[k], solRPdeIndex[k], i, iel); // global to global mapping between solution node and pdeSys dof } sysDof[3 * nDofs + i] = pdeSys->GetSystemDof(solHIndex, solHPdeIndex, i, iel); // global to global mapping between solution node and pdeSys dof } // local storage of coordinates for (unsigned i = 0; i < nDofs2; i++) { unsigned xDof = msh->GetSolutionDof(i, iel, xType); // global to global mapping between coordinates node and coordinate dof for (unsigned idim = 0; idim < dim; idim++) { x[idim][i] = (*msh->_topology->_Sol[idim])(xDof); // global extraction and local storage for the element coordinates } } // start a new recording of all the operations involving adept::adouble variables s.new_recording(); // *** Gauss point loop *** for (unsigned ig = 0; ig < msh->_finiteElement[ielGeom][solHType]->GetGaussPointNumber(); ig++) { // *** get gauss point weight, test function and test function partial derivatives *** msh->_finiteElement[ielGeom][solHType]->Jacobian(x, ig, weight, phi, phi_x, phi_xx); // evaluate the solution, the solution derivatives and the coordinates in the gauss point adept::adouble solRGauss[3]; double solRGaussOld[3]; adept::adouble solRGauss_x[3][2]; adept::adouble solRGauss_xx[3][2][2]; adept::adouble sol_x[2]; sol_x[0] = sol_x[1] = 0.; for (int k = 0; k < 3; k++) { solRGauss[k] = 0.; solRGaussOld[k] = 0.; for (int i = 0; i < dim; i++) { solRGauss_x[k][i] = 0.; for (int j = 0; j < dim; j++) { solRGauss_xx[k][i][j] = 0.; } } } adept::adouble solHGauss = 0; adept::adouble solHGauss_x[2] = {0., 0.}; for (unsigned i = 0; i < nDofs; i++) { for (int k = 0; k < 3; k++) { solRGauss[k] += phi[i] * solR[k][i]; solRGaussOld[k] += phi[i] * solR_old[k][i]; } solHGauss += phi[i] * solH[i]; for (unsigned u = 0; u < dim; u++) { sol_x[u] += phi[i] * x[u][i]; } for (unsigned u = 0; u < dim; u++) { // gradient for (int k = 0; k < 3; k++) { solRGauss_x[k][u] += phi_x[i * dim + u] * solR[k][i]; } solHGauss_x[u] += phi_x[i * dim + u] * solH[i]; } for (unsigned u = 0; u < dim; u++) { // hessian for (unsigned v = 0; v < dim; v++) { unsigned uvindex = 0; //_uu if (u != v) uvindex = 2; //_uv or _vu else if (u == 1) uvindex = 1; //_vv for (int k = 0; k < 3; k++) { solRGauss_xx[k][u][v] += phi_xx[i * dim2 + uvindex] * solR[k][i]; } } } } adept::adouble g[2][2]; g[0][0] = g[0][1] = g[1][0] = g[1][1] = 0.; for (int k = 0; k < 3; k++) { for (int u = 0; u < dim; u++) { for (int v = 0; v < dim; v++) { g[u][v] += solRGauss_x[k][u] * solRGauss_x[k][v]; } } } adept::adouble detg = g[0][0] * g[1][1] - g[0][1] * g[1][0]; adept::adouble A = sqrt(detg); adept::adouble gI[2][2]; gI[0][0] = g[1][1] / detg; gI[0][1] = -g[0][1] / detg; gI[1][0] = -g[1][0] / detg; gI[1][1] = g[0][0] / detg; adept::adouble N[3]; N[0] = (solRGauss_x[1][0] * solRGauss_x[2][1] - solRGauss_x[1][1] * solRGauss_x[2][0]) / A; N[1] = (solRGauss_x[2][0] * solRGauss_x[0][1] - solRGauss_x[2][1] * solRGauss_x[0][0]) / A; N[2] = (solRGauss_x[0][0] * solRGauss_x[1][1] - solRGauss_x[0][1] * solRGauss_x[1][0]) / A; adept::adouble h[2][2]; h[0][0] = h[0][1] = h[1][0] = h[1][1] = 0.; for (int k = 0; k < 3; k++) { for (int u = 0; u < dim; u++) { for (int v = 0; v < dim; v++) { h[u][v] += solRGauss_xx[k][u][v] * N[k]; } } } //adept::adouble K = cos(sol_x[0])/(a+cos(sol_x[0]));//(h[0][0]*h[1][1]-h[0][1]*h[1][0])/detg; adept::adouble K = (h[0][0] * h[1][1] - h[0][1] * h[1][0]) / detg; adept::adouble H_exact = 0.5 * (1. + cos(sol_x[0]) / (a + cos(sol_x[0]))); // *** phi_i loop *** for (unsigned i = 0; i < nDofs; i++) { for (int k = 0; k < 3; k++) { for (int u = 0; u < dim; u++) { adept::adouble AgIgradRgradPhi = 0; for (int v = 0; v < dim; v++) { AgIgradRgradPhi += A * gI[u][v].value() * solRGauss_x[k][v]; } aResR[k][i] += AgIgradRgradPhi * phi_x[i * dim + u] * weight; } aResR[k][i] += 2.* A * solHGauss.value() * N[k] * phi[i] * weight; } for (int u = 0; u < dim; u++) { adept::adouble AgIgradHgradPhi = 0; for (int v = 0; v < dim; v++) { AgIgradHgradPhi += A * gI[u][v].value() * solHGauss_x[v]; } aResH[i] -= AgIgradHgradPhi * phi_x[i * dim + u] * weight; } aResH[i] += A * (-0 * (solRGauss[0] - solRGaussOld[0]) * N[0].value() - 0 * (solRGauss[1] - solRGaussOld[1]) * N[1].value() - 0 * (solRGauss[2] - solRGaussOld[2]) * N[2].value() + 2. * solHGauss * (solHGauss * solHGauss - K.value())) * phi[i] * weight; } // end phi_i loop } // end gauss point loop //-------------------------------------------------------------------------------------------------------- // Add the local Matrix/Vector into the global Matrix/Vector //copy the value of the adept::adoube aRes in double Res and store for (int i = 0; i < nDofs; i++) { for (int k = 0; k < 3; k++) { Res[ k * nDofs + i] = -aResR[k][i].value(); } Res[ 3 * nDofs + i] = -aResH[i].value(); } RES->add_vector_blocked(Res, sysDof); Jac.resize((4 * nDofs) * (4 * nDofs)); // define the dependent variables for (int k = 0; k < 3; k++) { s.dependent(&aResR[k][0], nDofs); } s.dependent(&aResH[0], nDofs); // define the independent variables for (int k = 0; k < 3; k++) { s.independent(&solR[k][0], nDofs); } s.independent(&solH[0], nDofs); // get the jacobian matrix (ordered by row) s.jacobian(&Jac[0], true); KK->add_matrix_blocked(Jac, sysDof, sysDof); s.clear_independents(); s.clear_dependents(); } //end element loop for each process RES->close(); KK->close(); // ***************** END ASSEMBLY ******************* }
//------------------------------------------------------------------------------------------------------------ void AssembleMatrixRes_VC(MultiLevelProblem &ml_prob) { TransientNonlinearImplicitSystem* mlPdeSys = & ml_prob.get_system<TransientNonlinearImplicitSystem>("Timedep"); const unsigned level = mlPdeSys->GetLevelToAssemble(); MultiLevelSolution *ml_sol = ml_prob._ml_sol; Solution* sol = ml_sol->GetSolutionLevel(level); LinearEquationSolver* pdeSys = mlPdeSys->_LinSolver[level]; const char* pdename = mlPdeSys->name().c_str(); Mesh* msh = ml_prob._ml_msh->GetLevel(level); elem* myel = msh->el; SparseMatrix* JAC = pdeSys->_KK; NumericVector* RES = pdeSys->_RES; // data const unsigned dim = msh->GetDimension(); const unsigned max_size = static_cast< unsigned >(ceil(pow(3, dim))); // conservative: based on line3, quad9, hex27 unsigned nel = msh->GetNumberOfElements(); unsigned igrid = msh->GetLevel(); unsigned iproc = msh->processor_id(); // time dep data double dt = mlPdeSys->GetIntervalTime(); // double theta = 0.5; //************** geometry (at dofs and quadrature points) ************************************* vector < vector < double > > coords_at_dofs(dim); unsigned coords_fe_type = BIQUADR_FE; // get the finite element type for "x", it is always 2 (LAGRANGE BIQUADRATIC) for (unsigned i = 0; i < coords_at_dofs.size(); i++) coords_at_dofs[i].reserve(max_size); vector < double > coord_at_qp(dim); //************* shape functions (at dofs and quadrature points) ************************************** const int solType_max = BIQUADR_FE; //biquadratic double weight_qp; // gauss point weight vector < vector < double > > phi_fe_qp(NFE_FAMS); vector < vector < double > > phi_x_fe_qp(NFE_FAMS); vector < vector < double > > phi_xx_fe_qp(NFE_FAMS); for(int fe=0; fe < NFE_FAMS; fe++) { phi_fe_qp[fe].reserve(max_size); phi_x_fe_qp[fe].reserve(max_size*dim); phi_xx_fe_qp[fe].reserve(max_size*(3*(dim-1))); } //*************************************************** //********* WHOLE SET OF VARIABLES ****************** //*************************************************** const unsigned int n_unknowns = mlPdeSys->GetSolPdeIndex().size(); vector < std::string > Solname(n_unknowns); Solname[0] = "u"; vector < unsigned > SolPdeIndex(n_unknowns); vector < unsigned > SolIndex(n_unknowns); vector < unsigned int > SolFEType(n_unknowns); //FEtype of each MultilevelSolution vector < unsigned int > Sol_n_el_dofs(n_unknowns); //number of element dofs std::fill(Sol_n_el_dofs.begin(), Sol_n_el_dofs.end(), 0); for(unsigned ivar=0; ivar < n_unknowns; ivar++) { SolPdeIndex[ivar] = mlPdeSys->GetSolPdeIndex( Solname[ivar].c_str() ); SolIndex[ivar] = ml_sol->GetIndex ( Solname[ivar].c_str() ); SolFEType[ivar] = ml_sol->GetSolutionType ( SolIndex[ivar]); } //------------ quantities (at quadrature points) --------------------- vector<double> sol_qp(n_unknowns); vector<double> sol_old_qp(n_unknowns); vector< vector<double> > sol_grad_qp(n_unknowns); vector< vector<double> > sol_old_grad_qp(n_unknowns); std::fill(sol_qp.begin(), sol_qp.end(), 0.); std::fill(sol_old_qp.begin(), sol_old_qp.end(), 0.); for (unsigned k = 0; k < n_unknowns; k++) { sol_grad_qp[k].resize(dim); std::fill(sol_grad_qp[k].begin(), sol_grad_qp[k].end(), 0.); sol_old_grad_qp[k].resize(dim); std::fill(sol_old_grad_qp[k].begin(), sol_old_grad_qp[k].end(), 0.); } //----------- quantities (at dof objects) ------------------------------ vector < vector < double > > sol_eldofs(n_unknowns); vector < vector < double > > sol_old_eldofs(n_unknowns); for(int k=0; k<n_unknowns; k++) { sol_eldofs[k].reserve(max_size); sol_old_eldofs[k].reserve(max_size); } //******** linear system ******************************************* vector < vector < int > > L2G_dofmap(n_unknowns); for(int i = 0; i < n_unknowns; i++) { L2G_dofmap[i].reserve(max_size); } vector< int > L2G_dofmap_AllVars; L2G_dofmap_AllVars.reserve( n_unknowns*max_size ); vector< vector< double > > Res_el(n_unknowns); vector< vector< vector< double > > > Jac_el(n_unknowns); for(int i = 0; i < n_unknowns; i++) Res_el[i].reserve(max_size); for(int i = 0; i < n_unknowns; i++) { Jac_el[i].resize(n_unknowns); for(int j = 0; j < n_unknowns; j++) { Jac_el[i][j].reserve(max_size*max_size); } } // Set to zero all the entries of the matrix RES->zero(); JAC->zero(); const double deltat_term = 1.; const double lapl_term = 1.; const double delta_g_term = 1.; // *** element loop *** for (int iel = msh->_elementOffset[iproc]; iel < msh->_elementOffset[iproc+1]; iel++) { short unsigned ielGeom = msh->GetElementType(iel); // element geometry type //******************** GEOMETRY ********************* unsigned nDofx = msh->GetElementDofNumber(iel, coords_fe_type); // number of coordinate element dofs for (int i = 0; i < dim; i++) coords_at_dofs[i].resize(nDofx); // local storage of coordinates for (unsigned i = 0; i < nDofx; i++) { unsigned xDof = msh->GetSolutionDof(i, iel, coords_fe_type); // global to global mapping between coordinates node and coordinate dof for (unsigned jdim = 0; jdim < dim; jdim++) { coords_at_dofs[jdim][i] = (*msh->_topology->_Sol[jdim])(xDof); // global extraction and local storage for the element coordinates } } //*************************************************** //all vars################################################################### for (unsigned k = 0; k < n_unknowns; k++) { unsigned ndofs_unk = msh->GetElementDofNumber(iel, SolFEType[k]); Sol_n_el_dofs[k] = ndofs_unk; sol_eldofs[k].resize(ndofs_unk); sol_old_eldofs[k].resize(ndofs_unk); L2G_dofmap[k].resize(ndofs_unk); for (unsigned i = 0; i < ndofs_unk; i++) { unsigned solDof = msh->GetSolutionDof(i, iel, SolFEType[k]); // global to global mapping between solution node and solution dof // via local to global solution node sol_eldofs[k][i] = (*sol->_Sol[SolIndex[k]])(solDof); // global extraction and local storage for the solution sol_old_eldofs[k][i] = (*sol->_SolOld[SolIndex[k]])(solDof); // This is OLD in TIME, not in nonlinear loop L2G_dofmap[k][i] = pdeSys->GetSystemDof(SolIndex[k], SolPdeIndex[k], i, iel); // global to global mapping between solution node and pdeSys dof } } unsigned nDof_AllVars = 0; for (unsigned k = 0; k < n_unknowns; k++) { nDof_AllVars += Sol_n_el_dofs[k]; } // TODO COMPUTE MAXIMUM maximum number of element dofs for one scalar variable int nDof_max = 0; for (unsigned k = 0; k < n_unknowns; k++) { if(Sol_n_el_dofs[k] > nDof_max) nDof_max = Sol_n_el_dofs[k]; } for(int k = 0; k < n_unknowns; k++) { L2G_dofmap[k].resize(Sol_n_el_dofs[k]); Res_el[SolPdeIndex[k]].resize(Sol_n_el_dofs[k]); memset(& Res_el[SolPdeIndex[k]][0], 0., Sol_n_el_dofs[k] * sizeof(double) ); Jac_el[SolPdeIndex[k]][SolPdeIndex[k]].resize(Sol_n_el_dofs[k] * Sol_n_el_dofs[k]); memset(& Jac_el[SolPdeIndex[k]][SolPdeIndex[k]][0], 0., Sol_n_el_dofs[k] * Sol_n_el_dofs[k] * sizeof(double)); } //all vars################################################################### // *** quadrature loop *** for(unsigned ig = 0; ig < ml_prob.GetQuadratureRule(ielGeom).GetGaussPointsNumber(); ig++) { // *** get gauss point weight, test function and test function partial derivatives *** for(int fe=0; fe < NFE_FAMS; fe++) { msh->_finiteElement[ielGeom][fe]->Jacobian(coords_at_dofs,ig,weight_qp,phi_fe_qp[fe],phi_x_fe_qp[fe],phi_xx_fe_qp[fe]); } //HAVE TO RECALL IT TO HAVE BIQUADRATIC JACOBIAN msh->_finiteElement[ielGeom][coords_fe_type]->Jacobian(coords_at_dofs,ig,weight_qp,phi_fe_qp[coords_fe_type],phi_x_fe_qp[coords_fe_type],phi_xx_fe_qp[coords_fe_type]); //========= fill gauss value quantities ================== std::fill(sol_qp.begin(), sol_qp.end(), 0.); std::fill(sol_old_qp.begin(), sol_old_qp.end(), 0.); for (unsigned k = 0; k < n_unknowns; k++) { std::fill(sol_grad_qp[k].begin(), sol_grad_qp[k].end(), 0.); std::fill(sol_old_grad_qp[k].begin(), sol_old_grad_qp[k].end(), 0.); } for (unsigned k = 0; k < n_unknowns; k++) { for (unsigned i = 0; i < Sol_n_el_dofs[k]; i++) { sol_qp[k] += sol_eldofs[k][i] * phi_fe_qp[SolFEType[k]][i]; sol_old_qp[k] += sol_old_eldofs[k][i] * phi_fe_qp[SolFEType[k]][i]; for (unsigned d = 0; d < dim; d++) { sol_grad_qp[k][d] += sol_eldofs[k][i] * phi_x_fe_qp[SolFEType[k]][i * dim + d]; sol_old_grad_qp[k][d] += sol_old_eldofs[k][i] * phi_x_fe_qp[SolFEType[k]][i * dim + d]; } } } //========= fill gauss value quantities ================== // *** phi_i loop *** for(unsigned i = 0; i < nDof_max; i++) { //BEGIN RESIDUALS A block =========================== double Lap_rhs_i = 0.; double Lap_old_rhs_i = 0.; for(unsigned d = 0; d < dim; d++) { Lap_rhs_i += phi_x_fe_qp[SolFEType[0]][i * dim + d] * sol_grad_qp[0][d]; Lap_old_rhs_i += phi_x_fe_qp[SolFEType[0]][i * dim + d] * sol_old_grad_qp[0][d]; } Res_el[SolPdeIndex[0]][i] += - weight_qp * ( delta_g_term * Singularity::function(sol_qp[0]) * Singularity::g_vc(sol_qp[0]) * phi_fe_qp[ SolFEType[0] ][i] - delta_g_term * Singularity::function(sol_qp[0]) * Singularity::g_vc(sol_old_qp[0]) * phi_fe_qp[ SolFEType[0] ][i] + deltat_term * Singularity::function(sol_qp[0]) * dt * phi_fe_qp[ SolFEType[0] ][i] + lapl_term * dt * Lap_rhs_i ); //END RESIDUALS A block =========================== // *** phi_j loop *** for(unsigned j = 0; j<nDof_max; j++) { double Lap_mat_i_j = 0.; for(unsigned d = 0; d < dim; d++) Lap_mat_i_j += phi_x_fe_qp[SolFEType[0]][i * dim + d] * phi_x_fe_qp[SolFEType[0]][j * dim + d]; Jac_el[SolPdeIndex[0]][SolPdeIndex[0]][i * Sol_n_el_dofs[0] + j] += weight_qp * ( + delta_g_term * phi_fe_qp[ SolFEType[0] ][i] * Singularity::derivative( phi_fe_qp[ SolFEType[0] ][j] ) * phi_fe_qp[ SolFEType[0] ][j] * Singularity::g_vc( phi_fe_qp[ SolFEType[0] ][j] ) + delta_g_term * phi_fe_qp[ SolFEType[0] ][i] * Singularity::function ( phi_fe_qp[ SolFEType[0] ][j] ) * Singularity::g_vc_derivative( phi_fe_qp[ SolFEType[0] ][j] ) * phi_fe_qp[ SolFEType[0] ][j] - delta_g_term * phi_fe_qp[ SolFEType[0] ][i] * Singularity::derivative( phi_fe_qp[ SolFEType[0] ][j] ) * phi_fe_qp[ SolFEType[0] ][j] * Singularity::g_vc( sol_old_qp[0] ) + deltat_term * phi_fe_qp[ SolFEType[0] ][i] * Singularity::derivative( phi_fe_qp[ SolFEType[0] ][j] ) * phi_fe_qp[ SolFEType[0] ][j] * dt + lapl_term * dt * Lap_mat_i_j ); } //end phij loop } //end phii loop } // end gauss point loop //-------------------------------------------------------------------------------------------------------- //Sum the local matrices/vectors into the Global Matrix/Vector for(unsigned ivar=0; ivar < n_unknowns; ivar++) { RES->add_vector_blocked(Res_el[SolPdeIndex[ivar]],L2G_dofmap[ivar]); JAC->add_matrix_blocked(Jac_el[SolPdeIndex[ivar]][SolPdeIndex[ivar]],L2G_dofmap[ivar],L2G_dofmap[ivar]); } //-------------------------------------------------------------------------------------------------------- } //end list of elements loop for each subdomain JAC->close(); RES->close(); // ***************** END ASSEMBLY ******************* }