/** * This function assemble the stiffnes matrix KK and the residual vector Res * Using automatic differentiation for Newton iterative scheme * J(u0) w = - F(u0) , * with u = u0 + w * - F = f(x) - J u = Res * J = \grad_u F * * thus * J w = f(x) - J u0 **/ void AssemblePoissonProblem_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 LinearImplicitSystem* mlPdeSys = &ml_prob.get_system<LinearImplicitSystem> ("Poisson"); // pointer to the linear implicit system named "Poisson" const unsigned level = mlPdeSys->GetLevelToAssemble(); const unsigned levelMax = mlPdeSys->GetLevelMax(); const bool assembleMatrix = mlPdeSys->GetAssembleMatrix(); 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 dim2 = (3 * (dim - 1) + !(dim - 1)); // dim2 is the number of second order partial derivatives (1,3,6 depending on the dimension) const unsigned maxSize = static_cast< unsigned >(ceil(pow(3, dim))); // conservative: based on line3, quad9, hex27 unsigned iproc = msh->processor_id(); // get the process_id (for parallel computation) //solution variable unsigned soluIndex; soluIndex = mlSol->GetIndex("u"); // get the position of "u" in the ml_sol object unsigned soluType = mlSol->GetSolutionType(soluIndex); // get the finite element type for "u" unsigned soluPdeIndex; soluPdeIndex = mlPdeSys->GetSolPdeIndex("u"); // get the position of "u" in the pdeSys object vector < adept::adouble > solu; // local solution solu.reserve(maxSize); vector < vector < double > > x(dim); // local coordinates unsigned xType = 2; // get the finite element type for "x", it is always 2 (LAGRANGE QUADRATIC) for (unsigned i = 0; i < dim; i++) { x[i].reserve(maxSize); } 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 phi.reserve(maxSize); phi_x.reserve(maxSize * dim); phi_xx.reserve(maxSize * dim2); vector< adept::adouble > aRes; // local redidual vector aRes.reserve(maxSize); vector< int > l2GMap; // local to global mapping l2GMap.reserve(maxSize); vector< double > Res; // local redidual vector Res.reserve(maxSize); vector < double > Jac; Jac.reserve(maxSize * 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->IS_Mts2Gmt_elem_offset[iproc]; iel < msh->IS_Mts2Gmt_elem_offset[iproc + 1]; iel++) { unsigned kel = msh->IS_Mts2Gmt_elem[iel]; // mapping between paralell dof and mesh dof short unsigned kelGeom = el->GetElementType(kel); // element geometry type unsigned nDofu = el->GetElementDofNumber(kel, soluType); // number of solution element dofs unsigned nDofx = el->GetElementDofNumber(kel, xType); // number of coordinate element dofs // resize local arrays l2GMap.resize(nDofu); solu.resize(nDofu); for (int i = 0; i < dim; i++) { x[i].resize(nDofx); } aRes.resize(nDofu); //resize std::fill(aRes.begin(), aRes.end(), 0); //set aRes to zero // local storage of global mapping and solution for (unsigned i = 0; i < nDofu; i++) { unsigned iNode = el->GetMeshDof(kel, i, soluType); // local to global solution node unsigned solDof = msh->GetMetisDof(iNode, soluType); // global to global mapping between solution node and solution dof solu[i] = (*sol->_Sol[soluIndex])(solDof); // global extraction and local storage for the solution l2GMap[i] = pdeSys->GetKKDof(soluIndex, soluPdeIndex, iNode); // global to global mapping between solution node and pdeSys dof } // local storage of coordinates for (unsigned i = 0; i < nDofx; i++) { unsigned iNode = el->GetMeshDof(kel, i, xType); // local to global coordinates node unsigned xDof = msh->GetMetisDof(iNode, xType); // global to global mapping between coordinates node and coordinate dof for (unsigned jdim = 0; jdim < dim; jdim++) { x[jdim][i] = (*msh->_coordinate->_Sol[jdim])(xDof); // global extraction and local storage for the element coordinates } } if (level == levelMax || !el->GetRefinedElementIndex(kel)) { // do not care about this if now (it is used for the AMR) // 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[kelGeom][soluType]->GetGaussPointNumber(); ig++) { // *** get gauss point weight, test function and test function partial derivatives *** msh->_finiteElement[kelGeom][soluType]->Jacobian(x, ig, weight, phi, phi_x, phi_xx); // evaluate the solution, the solution derivatives and the coordinates in the gauss point adept::adouble solu_gss = 0; vector < adept::adouble > gradSolu_gss(dim, 0.); vector < double > x_gss(dim, 0.); for (unsigned i = 0; i < nDofu; i++) { solu_gss += phi[i] * solu[i]; for (unsigned jdim = 0; jdim < dim; jdim++) { gradSolu_gss[jdim] += phi_x[i * dim + jdim] * solu[i]; x_gss[jdim] += x[jdim][i] * phi[i]; } } // *** phi_i loop *** for (unsigned i = 0; i < nDofu; i++) { adept::adouble laplace = 0.; for (unsigned jdim = 0; jdim < dim; jdim++) { laplace += phi_x[i * dim + jdim] * gradSolu_gss[jdim]; } double srcTerm = - GetExactSolutionLaplace(x_gss); aRes[i] += (srcTerm * phi[i] - laplace) * weight; } // end phi_i loop } // end gauss point loop } // endif single element not refined or fine grid loop //-------------------------------------------------------------------------------------------------------- // Add the local Matrix/Vector into the global Matrix/Vector //copy the value of the adept::adoube aRes in double Res and store Res.resize(nDofu); //resize for (int i = 0; i < nDofu; i++) { Res[i] = - aRes[i].value(); } RES->add_vector_blocked(Res, l2GMap); if (assembleMatrix) { // define the dependent variables s.dependent(&aRes[0], nDofu); // define the independent variables s.independent(&solu[0], nDofu); // get the jacobian matrix (ordered by row major ) Jac.resize(nDofu * nDofu); //resize s.jacobian(&Jac[0], true); //store K in the global matrix KK KK->add_matrix_blocked(Jac, l2GMap, l2GMap); s.clear_independents(); s.clear_dependents(); } } //end element loop for each process RES->close(); if (assembleMatrix) KK->close(); // ***************** END ASSEMBLY ******************* }
void ETD(MultiLevelProblem& ml_prob) { const unsigned& NLayers = NumberOfLayers; adept::Stack& s = FemusInit::_adeptStack; LinearImplicitSystem* mlPdeSys = &ml_prob.get_system<LinearImplicitSystem> ("SW"); // pointer to the linear implicit system named "Poisson" unsigned level = ml_prob._ml_msh->GetNumberOfLevels() - 1u; 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) NumericVector* EPS = pdeSys->_EPS; // 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) unsigned nprocs = msh->n_processors(); // get the process_id (for parallel computation) //solution variable std::vector < unsigned > solIndexh(NLayers); std::vector < unsigned > solPdeIndexh(NLayers); std::vector < unsigned > solIndexv(NLayers); std::vector < unsigned > solPdeIndexv(NLayers); // std::vector < unsigned > solIndextracer(NLayers); // std::vector < unsigned > solPdeIndextracer(NLayers); vector< int > l2GMap; // local to global mapping for(unsigned i = 0; i < NLayers; i++) { char name[10]; sprintf(name, "h%d", i); solIndexh[i] = mlSol->GetIndex(name); // get the position of "hi" in the sol object solPdeIndexh[i] = mlPdeSys->GetSolPdeIndex(name); // get the position of "hi" in the pdeSys object sprintf(name, "v%d", i); solIndexv[i] = mlSol->GetIndex(name); // get the position of "vi" in the sol object solPdeIndexv[i] = mlPdeSys->GetSolPdeIndex(name); // get the position of "vi" in the pdeSys object // sprintf(name, "tracer%d", i); // solIndextracer[i] = mlSol->GetIndex(name); // get the position of "traceri" in the sol object // solPdeIndextracer[i] = mlPdeSys->GetSolPdeIndex(name); // get the position of "traceri" in the pdeSys object } unsigned solTypeh = mlSol->GetSolutionType(solIndexh[0]); // get the finite element type for "hi" unsigned solTypev = mlSol->GetSolutionType(solIndexv[0]); // get the finite element type for "vi" // unsigned solTypetracer = mlSol->GetSolutionType(solIndextracer[0]); // get the finite element type for "traceri" vector < double > x; // local coordinates vector< vector < adept::adouble > > solh(NLayers); // local coordinates vector< vector < adept::adouble > > solv(NLayers); // local coordinates // vector< vector < adept::adouble > > soltracer(NLayers); // local coordinates vector< vector < bool > > bdch(NLayers); // local coordinates vector< vector < bool > > bdcv(NLayers); // local coordinates // vector< vector < bool > > bdctracer(NLayers); // local coordinates unsigned xType = 2; // get the finite element type for "x", it is always 2 (LAGRANGE QUADRATIC) vector < vector< adept::adouble > > aResh(NLayers); vector < vector< adept::adouble > > aResv(NLayers); // vector < vector< adept::adouble > > aRestracer(NLayers); KK->zero(); RES->zero(); double maxWaveSpeed = 0.; double dx; for(int iel = msh->_elementOffset[iproc]; iel < msh->_elementOffset[iproc + 1]; iel++) { short unsigned ielGeom = msh->GetElementType(iel); unsigned nDofh = msh->GetElementDofNumber(iel, solTypeh); // number of solution element dofs unsigned nDofv = msh->GetElementDofNumber(iel, solTypev); // number of solution element dofs //unsigned nDoftracer = msh->GetElementDofNumber(iel, solTypetracer); // number of coordinate element dofs unsigned nDofx = msh->GetElementDofNumber(iel, xType); // number of coordinate element dofs unsigned nDofs = nDofh + nDofv /*+ nDoftracer*/; // resize local arrays l2GMap.resize(NLayers * (nDofs) ); for(unsigned i = 0; i < NLayers; i++) { solh[i].resize(nDofh); solv[i].resize(nDofv); //soltracer[i].resize(nDofvtracer); bdch[i].resize(nDofh); bdcv[i].resize(nDofv); //bdctracer[i].resize(nDoftracer); aResh[i].resize(nDofh); //resize std::fill(aResh[i].begin(), aResh[i].end(), 0); //set aRes to zero aResv[i].resize(nDofv); //resize std::fill(aResv[i].begin(), aResv[i].end(), 0); //set aRes to zero //aRestracer[i].resize(nDoftracer); //resize //std::fill(aRestracer[i].begin(), aRestracer[i].end(), 0); //set aRes to zero } x.resize(nDofx); //local storage of global mapping and solution for(unsigned i = 0; i < nDofh; i++) { unsigned solDofh = msh->GetSolutionDof(i, iel, solTypeh); // global to global mapping between solution node and solution dof for(unsigned j = 0; j < NLayers; j++) { solh[j][i] = (*sol->_Sol[solIndexh[j]])(solDofh); // global extraction and local storage for the solution bdch[j][i] = ( (*sol->_Bdc[solIndexh[j]])(solDofh) < 1.5) ? true : false; l2GMap[ j * nDofs + i] = pdeSys->GetSystemDof(solIndexh[j], solPdeIndexh[j], i, iel); // global to global mapping between solution node and pdeSys dof } } for(unsigned i = 0; i < nDofv; i++) { unsigned solDofv = msh->GetSolutionDof(i, iel, solTypev); // global to global mapping between solution node and solution dof for(unsigned j = 0; j < NLayers; j++) { solv[j][i] = (*sol->_Sol[solIndexv[j]])(solDofv); // global extraction and local storage for the solution bdcv[j][i] = ( (*sol->_Bdc[solIndexv[j]])(solDofv) < 1.5) ? true : false; l2GMap[ j * nDofs + nDofh + i] = pdeSys->GetSystemDof(solIndexv[j], solPdeIndexv[j], i, iel); // global to global mapping between solution node and pdeSys dof } } // for(unsigned i = 0; i < nDoftracer; i++) { // unsigned solDoftracer = msh->GetSolutionDof(i, iel, solTypetracer); // global to global mapping between solution node and solution dof // for(unsigned j = 0; j < NLayers; j++) { // soltracer[j][i] = (*sol->_Sol[solIndextracer[j]])(solDoftracer); // global extraction and local storage for the solution // bdctracer[j][i] = ( (*sol->_Bdc[solIndextracer[j]])(solDoftracer) < 1.5) ? true : false; // l2GMap[ j * nDofs + nDofh + nDofv + i] = pdeSys->GetSystemDof(solIndextracer[j], solPdeIndextracer[j], i, iel); // global to global mapping between solution node and pdeSys dof // } // } s.new_recording(); for(unsigned i = 0; i < nDofx; i++) { unsigned xDof = msh->GetSolutionDof(i, iel, xType); // global to global mapping between coordinates node and coordinate dof x[i] = (*msh->_topology->_Sol[0])(xDof); // global extraction and local storage for the element coordinates } double xmid = 0.5 * (x[0] + x[1]); //std::cout << xmid << std::endl; double zz = sqrt(aa * aa - xmid * xmid); // z coordinate of points on sphere double dd = aa * acos((zz * z_c) / (aa * aa)); // distance to center point on sphere [m] double hh = 1 - dd * dd / (bb * bb); double b = ( H_shelf + H_0 / 2 * (1 + tanh(hh / phi)) ); std::vector < adept::adouble > eta(NLayers); double hTot = 0.; for(unsigned k = 0; k < NLayers; k++) { hTot += solh[k][0].value(); eta[k] = - b * rho[k]; // bottom topography for( unsigned j = 0; j < NLayers; j++){ double rhoj = (j <= k) ? rho[j] : rho[k]; eta[k] += rhoj * solh[j][0]; } eta[k] /= rho[k]; } std::vector < double > hALE(NLayers, 0.); hALE[0] = hRest[0] + (hTot - b); for(unsigned k = 1; k < NLayers - 1; k++){ hALE[k] = hRest[k]; } hALE[NLayers - 1] = b - hRest[NLayers - 1]; std::vector < double > w(NLayers+1, 0.); dx = x[1] - x[0]; for(unsigned k = NLayers; k>1; k--){ w[k-1] = w[k] - solh[k-1][0].value() * (solv[k-1][1].value() - solv[k-1][0].value() )/dx- ( hALE[k-1] - solh[k-1][0].value()) / dt; //std::cout << hALE[k-1] << " " << w[k-1] << " "; } //std::cout<<std::endl; std::vector < adept::adouble > zMid(NLayers); for(unsigned k = 0; k < NLayers; k++) { zMid[k] = solh[k][0]/2; for(unsigned i = 0; i < k; i++) { zMid[k] += solh[i][0]; } //std::cout << zMid[k] << " "; } //std::cout<<std::endl; // double hTotal = 0.; // for(unsigned i = 0; i<NLayers; i++) hTotal += solh[i][0].value(); double celerity = sqrt( 9.81 * hTot); for(unsigned i = 0; i<NLayers; i++){ double vmid = 0.5 * ( solv[i][0].value() + solv[i][1].value() ); maxWaveSpeed = ( maxWaveSpeed > fabs(vmid) + fabs(celerity) )? maxWaveSpeed : fabs(vmid) + fabs(celerity); } for(unsigned k = 0; k < NLayers; k++) { if(!bdch[k][0]) { for (unsigned j = 0; j < nDofv; j++) { double sign = ( j == 0) ? 1. : -1; aResh[k][0] += sign * solh[k][0] * solv[k][j] / dx; } aResh[k][0] += w[k+1] - w[k]; } adept::adouble vMid = 0.5 * (solv[k][0] + solv[k][1]); adept::adouble fv = 0.5 * vMid * vMid + 9.81 * eta[k]; for (unsigned i = 0; i < nDofv; i++) { if(!bdcv[k][i]) { double sign = ( i == 0) ? -1. : 1; aResv[k][i] += sign * fv / dx; if ( k > 0 ){ aResv[k][i] -= 0.5 * ( w[k] * 0.5 * ( solv[k-1][i] - solv[k][i] ) / (0.5 * ( solh[k-1][0] + solh[k][0] ) ) ); } if (k < NLayers - 1) { aResv[k][i] -= 0.5 * ( w[k+1] * 0.5 * ( solv[k][i] - solv[k+1][i] ) / (0.5 * ( solh[k][0] + solh[k+1][0] ) ) ); } //aResv[k][i] += sign * 9.81 * rho[k] * zMid[k] / dx; } } // if(!bdctracer[k][0]) { // for (unsigned j = 0; j < nDofv; j++) { // double sign = ( j == 0) ? 1. : -1; // aRestracer[k][0] += sign * solh[k][0] * solv[k][j] * soltracer[k][0] / dx; // } // aRestracer[k][0] += w[k+1] * ((soltracer[k][0] + soltracer[k+1][0])/2) - w[k] * ((soltracer[k-1][0] + soltracer[k][0])/2); // } } vector< double > Res(NLayers * nDofs); // local redidual vector unsigned counter = 0; for(unsigned k = 0; k < NLayers; k++) { for(int i = 0; i < nDofh; i++) { Res[counter] = aResh[k][i].value(); counter++; } for(int i = 0; i < nDofv; i++) { Res[counter] = aResv[k][i].value(); counter++; } } RES->add_vector_blocked(Res, l2GMap); for(unsigned k = 0; k < NLayers; k++) { // define the dependent variables s.dependent(&aResh[k][0], nDofh); s.dependent(&aResv[k][0], nDofv); // define the independent variables s.independent(&solh[k][0], nDofh); s.independent(&solv[k][0], nDofv); } // get the jacobian matrix (ordered by row major ) vector < double > Jac(NLayers * nDofs * NLayers * nDofs); s.jacobian(&Jac[0], true); //store K in the global matrix KK KK->add_matrix_blocked(Jac, l2GMap, l2GMap); s.clear_independents(); s.clear_dependents(); } RES->close(); KK->close(); // PetscViewer viewer; // PetscViewerDrawOpen(PETSC_COMM_WORLD,NULL,NULL,0,0,900,900,&viewer); // PetscObjectSetName((PetscObject)viewer,"FSI matrix"); // PetscViewerPushFormat(viewer,PETSC_VIEWER_DRAW_LG); // MatView((static_cast<PetscMatrix*>(KK))->mat(),viewer); // double a; // std::cin>>a; // // abort(); MFN mfn; Mat A = (static_cast<PetscMatrix*>(KK))->mat(); FN f, f1, f2, f3 , f4; std::cout << "dt = " << dt << " dx = "<< dx << " maxWaveSpeed = "<<maxWaveSpeed << std::endl; //dt = 100.; Vec v = (static_cast< PetscVector* >(RES))->vec(); Vec y = (static_cast< PetscVector* >(EPS))->vec(); MFNCreate( PETSC_COMM_WORLD, &mfn ); MFNSetOperator( mfn, A ); MFNGetFN( mfn, &f ); // FNCreate(PETSC_COMM_WORLD, &f1); // FNCreate(PETSC_COMM_WORLD, &f2); // FNCreate(PETSC_COMM_WORLD, &f3); // FNCreate(PETSC_COMM_WORLD, &f4); // // FNSetType(f1, FNEXP); // // FNSetType(f2, FNRATIONAL); // double coeff1[1] = { -1}; // FNRationalSetNumerator(f2, 1, coeff1); // FNRationalSetDenominator(f2, 0, PETSC_NULL); // // FNSetType( f3, FNCOMBINE ); // // FNCombineSetChildren(f3, FN_COMBINE_ADD, f1, f2); // // FNSetType(f4, FNRATIONAL); // double coeff2[2] = {1., 0.}; // FNRationalSetNumerator(f4, 2, coeff2); // FNRationalSetDenominator(f4, 0, PETSC_NULL); // // FNSetType( f, FNCOMBINE ); // // FNCombineSetChildren(f, FN_COMBINE_DIVIDE, f3, f4); FNPhiSetIndex(f,1); FNSetType( f, FNPHI ); // FNView(f,PETSC_VIEWER_STDOUT_WORLD); FNSetScale( f, dt, dt); MFNSetFromOptions( mfn ); MFNSolve( mfn, v, y); MFNDestroy( &mfn ); // FNDestroy(&f1); // FNDestroy(&f2); // FNDestroy(&f3); // FNDestroy(&f4); sol->UpdateSol(mlPdeSys->GetSolPdeIndex(), EPS, pdeSys->KKoffset); unsigned solIndexeta = mlSol->GetIndex("eta"); unsigned solIndexb = mlSol->GetIndex("b"); sol->_Sol[solIndexeta]->zero(); for(unsigned k=0;k<NumberOfLayers;k++){ sol->_Sol[solIndexeta]->add(*sol->_Sol[solIndexh[k]]); } sol->_Sol[solIndexeta]->add(-1,*sol->_Sol[solIndexb]); }