13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1177841947SLeila Ghaffari #include "../navierstokes.h" 1277841947SLeila Ghaffari #include "../qfunctions/mass.h" 1377841947SLeila Ghaffari 1477841947SLeila Ghaffari // Compute mass matrix for explicit scheme 1577841947SLeila Ghaffari PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, 1677841947SLeila Ghaffari Vec M) { 1777841947SLeila Ghaffari Vec M_loc; 1877841947SLeila Ghaffari CeedQFunction qf_mass; 1977841947SLeila Ghaffari CeedOperator op_mass; 2077841947SLeila Ghaffari CeedVector m_ceed, ones_vec; 2177841947SLeila Ghaffari CeedInt num_comp_q, q_data_size; 2277841947SLeila Ghaffari PetscErrorCode ierr; 2377841947SLeila Ghaffari PetscFunctionBeginUser; 2477841947SLeila Ghaffari 2577841947SLeila Ghaffari // CEED Restriction 2677841947SLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q); 2777841947SLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size); 2877841947SLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL); 2977841947SLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL); 3077841947SLeila Ghaffari CeedVectorSetValue(ones_vec, 1.0); 3177841947SLeila Ghaffari 3277841947SLeila Ghaffari // CEED QFunction 3377841947SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 3477841947SLeila Ghaffari CeedQFunctionAddInput(qf_mass, "q", num_comp_q, CEED_EVAL_INTERP); 35a61c78d6SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE); 3677841947SLeila Ghaffari CeedQFunctionAddOutput(qf_mass, "v", num_comp_q, CEED_EVAL_INTERP); 3777841947SLeila Ghaffari 3877841947SLeila Ghaffari // CEED Operator 3977841947SLeila Ghaffari CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 4077841947SLeila Ghaffari CeedOperatorSetField(op_mass, "q", ceed_data->elem_restr_q, ceed_data->basis_q, 4177841947SLeila Ghaffari CEED_VECTOR_ACTIVE); 42a61c78d6SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, 4377841947SLeila Ghaffari CEED_BASIS_COLLOCATED, ceed_data->q_data); 4477841947SLeila Ghaffari CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 4577841947SLeila Ghaffari CEED_VECTOR_ACTIVE); 4677841947SLeila Ghaffari 4777841947SLeila Ghaffari // Place PETSc vector in CEED vector 4877841947SLeila Ghaffari CeedScalar *m; 4977841947SLeila Ghaffari PetscMemType m_mem_type; 5077841947SLeila Ghaffari ierr = DMGetLocalVector(dm, &M_loc); CHKERRQ(ierr); 5177841947SLeila Ghaffari ierr = VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type); 5277841947SLeila Ghaffari CHKERRQ(ierr); 5377841947SLeila Ghaffari CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m); 5477841947SLeila Ghaffari 5577841947SLeila Ghaffari // Apply CEED Operator 5677841947SLeila Ghaffari CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE); 5777841947SLeila Ghaffari 5877841947SLeila Ghaffari // Restore vectors 5977841947SLeila Ghaffari CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL); 6077841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m); 6177841947SLeila Ghaffari CHKERRQ(ierr); 6277841947SLeila Ghaffari 6377841947SLeila Ghaffari // Local-to-Global 6477841947SLeila Ghaffari ierr = VecZeroEntries(M); CHKERRQ(ierr); 6577841947SLeila Ghaffari ierr = DMLocalToGlobal(dm, M_loc, ADD_VALUES, M); CHKERRQ(ierr); 6677841947SLeila Ghaffari ierr = DMRestoreLocalVector(dm, &M_loc); CHKERRQ(ierr); 6777841947SLeila Ghaffari 6877841947SLeila Ghaffari // Invert diagonally lumped mass vector for RHS function 6977841947SLeila Ghaffari ierr = VecReciprocal(M); CHKERRQ(ierr); 7077841947SLeila Ghaffari 7177841947SLeila Ghaffari // Cleanup 7277841947SLeila Ghaffari CeedVectorDestroy(&ones_vec); 7377841947SLeila Ghaffari CeedVectorDestroy(&m_ceed); 7477841947SLeila Ghaffari CeedQFunctionDestroy(&qf_mass); 7577841947SLeila Ghaffari CeedOperatorDestroy(&op_mass); 7677841947SLeila Ghaffari 7777841947SLeila Ghaffari PetscFunctionReturn(0); 7877841947SLeila Ghaffari } 7977841947SLeila Ghaffari 8077841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup 8177841947SLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 8277841947SLeila Ghaffari // This function takes in a state vector Q and writes into G 8377841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 8477841947SLeila Ghaffari User user = *(User *)user_data; 8577841947SLeila Ghaffari PetscScalar *q, *g; 865e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_loc, G_loc; 8777841947SLeila Ghaffari PetscMemType q_mem_type, g_mem_type; 8877841947SLeila Ghaffari PetscErrorCode ierr; 8977841947SLeila Ghaffari PetscFunctionBeginUser; 9077841947SLeila Ghaffari 915e82a6e1SJeremy L Thompson // Get local vector 925e82a6e1SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 935e82a6e1SJeremy L Thompson 945e82a6e1SJeremy L Thompson // Update time dependent data 955e82a6e1SJeremy L Thompson if (user->time != t) { 965e82a6e1SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, 975e82a6e1SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 985e82a6e1SJeremy L Thompson if (user->phys->solution_time_label) { 9911436a05SJames Wright CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t); 1005e82a6e1SJeremy L Thompson } 1015e82a6e1SJeremy L Thompson user->time = t; 1025e82a6e1SJeremy L Thompson } 10388626eedSJames Wright if (user->phys->timestep_size_label) { 10488626eedSJames Wright PetscScalar dt; 10588626eedSJames Wright ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr); 1065e82a6e1SJeremy L Thompson if (user->dt != dt) { 10788626eedSJames Wright CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label, 10888626eedSJames Wright &dt); 1095e82a6e1SJeremy L Thompson user->dt = dt; 1105e82a6e1SJeremy L Thompson } 11188626eedSJames Wright } 11277841947SLeila Ghaffari 11377841947SLeila Ghaffari // Global-to-local 11477841947SLeila Ghaffari ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr); 11577841947SLeila Ghaffari 11677841947SLeila Ghaffari // Place PETSc vectors in CEED vectors 11777841947SLeila Ghaffari ierr = VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type); 11877841947SLeila Ghaffari CHKERRQ(ierr); 11977841947SLeila Ghaffari ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr); 12077841947SLeila Ghaffari CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q); 12177841947SLeila Ghaffari CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 12277841947SLeila Ghaffari 12377841947SLeila Ghaffari // Apply CEED operator 12477841947SLeila Ghaffari CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, 12577841947SLeila Ghaffari CEED_REQUEST_IMMEDIATE); 12677841947SLeila Ghaffari 12777841947SLeila Ghaffari // Restore vectors 12877841947SLeila Ghaffari CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 12977841947SLeila Ghaffari CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 13077841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q); 13177841947SLeila Ghaffari CHKERRQ(ierr); 13277841947SLeila Ghaffari ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr); 13377841947SLeila Ghaffari 13477841947SLeila Ghaffari // Local-to-Global 13577841947SLeila Ghaffari ierr = VecZeroEntries(G); CHKERRQ(ierr); 13677841947SLeila Ghaffari ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr); 13777841947SLeila Ghaffari 13877841947SLeila Ghaffari // Inverse of the lumped mass matrix (M is Minv) 13977841947SLeila Ghaffari ierr = VecPointwiseMult(G, G, user->M); CHKERRQ(ierr); 14077841947SLeila Ghaffari 14177841947SLeila Ghaffari // Restore vectors 14277841947SLeila Ghaffari ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 14377841947SLeila Ghaffari 14477841947SLeila Ghaffari PetscFunctionReturn(0); 14577841947SLeila Ghaffari } 14677841947SLeila Ghaffari 14777841947SLeila Ghaffari // Implicit time-stepper function setup 14877841947SLeila Ghaffari PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, 14977841947SLeila Ghaffari void *user_data) { 15077841947SLeila Ghaffari User user = *(User *)user_data; 15177841947SLeila Ghaffari const PetscScalar *q, *q_dot; 15277841947SLeila Ghaffari PetscScalar *g; 1535e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 15477841947SLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 15577841947SLeila Ghaffari PetscErrorCode ierr; 15677841947SLeila Ghaffari PetscFunctionBeginUser; 15777841947SLeila Ghaffari 1585e82a6e1SJeremy L Thompson // Get local vectors 1595e82a6e1SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 1605e82a6e1SJeremy L Thompson 1615e82a6e1SJeremy L Thompson // Update time dependent data 1625e82a6e1SJeremy L Thompson if (user->time != t) { 1635e82a6e1SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, 1645e82a6e1SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 1655e82a6e1SJeremy L Thompson if (user->phys->solution_time_label) { 16611436a05SJames Wright CeedOperatorContextSetDouble(user->op_ifunction, 16711436a05SJames Wright user->phys->solution_time_label, &t); 1685e82a6e1SJeremy L Thompson } 1695e82a6e1SJeremy L Thompson user->time = t; 1705e82a6e1SJeremy L Thompson } 17188626eedSJames Wright if (user->phys->timestep_size_label) { 17288626eedSJames Wright PetscScalar dt; 17388626eedSJames Wright ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr); 1745e82a6e1SJeremy L Thompson if (user->dt != dt) { 17588626eedSJames Wright CeedOperatorContextSetDouble(user->op_ifunction, 17688626eedSJames Wright user->phys->timestep_size_label, &dt); 1775e82a6e1SJeremy L Thompson user->dt = dt; 1785e82a6e1SJeremy L Thompson } 17988626eedSJames Wright } 18077841947SLeila Ghaffari 18177841947SLeila Ghaffari // Global-to-local 18277841947SLeila Ghaffari ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr); 18377841947SLeila Ghaffari ierr = DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc); 18477841947SLeila Ghaffari CHKERRQ(ierr); 18577841947SLeila Ghaffari 18677841947SLeila Ghaffari // Place PETSc vectors in CEED vectors 18777841947SLeila Ghaffari ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr); 18877841947SLeila Ghaffari ierr = VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type); 18977841947SLeila Ghaffari CHKERRQ(ierr); 19077841947SLeila Ghaffari ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr); 19177841947SLeila Ghaffari CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, 19277841947SLeila Ghaffari (PetscScalar *)q); 19377841947SLeila Ghaffari CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), 194e6225c47SLeila Ghaffari CEED_USE_POINTER, (PetscScalar *)q_dot); 19577841947SLeila Ghaffari CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 19677841947SLeila Ghaffari 19777841947SLeila Ghaffari // Apply CEED operator 19877841947SLeila Ghaffari CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, 19977841947SLeila Ghaffari CEED_REQUEST_IMMEDIATE); 20077841947SLeila Ghaffari 20177841947SLeila Ghaffari // Restore vectors 20277841947SLeila Ghaffari CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 20377841947SLeila Ghaffari CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL); 20477841947SLeila Ghaffari CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 20577841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr); 20677841947SLeila Ghaffari ierr = VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot); CHKERRQ(ierr); 20777841947SLeila Ghaffari ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr); 20877841947SLeila Ghaffari 20977841947SLeila Ghaffari // Local-to-Global 21077841947SLeila Ghaffari ierr = VecZeroEntries(G); CHKERRQ(ierr); 21177841947SLeila Ghaffari ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr); 21277841947SLeila Ghaffari 21377841947SLeila Ghaffari // Restore vectors 21477841947SLeila Ghaffari ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 21577841947SLeila Ghaffari 21677841947SLeila Ghaffari PetscFunctionReturn(0); 21777841947SLeila Ghaffari } 21877841947SLeila Ghaffari 219e334ad8fSJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) { 220e334ad8fSJed Brown User user; 221e334ad8fSJed Brown const PetscScalar *q; 222e334ad8fSJed Brown PetscScalar *g; 223e334ad8fSJed Brown PetscMemType q_mem_type, g_mem_type; 224e334ad8fSJed Brown PetscErrorCode ierr; 225e334ad8fSJed Brown PetscFunctionBeginUser; 2265e82a6e1SJeremy L Thompson ierr = MatShellGetContext(J, &user); CHKERRQ(ierr); 2275e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_dot_loc, // Note - Q_dot_loc has zero BCs 2285e82a6e1SJeremy L Thompson G_loc; 2295e82a6e1SJeremy L Thompson 230e334ad8fSJed Brown // Get local vectors 231e334ad8fSJed Brown ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 232e334ad8fSJed Brown 233e334ad8fSJed Brown // Global-to-local 234e334ad8fSJed Brown ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr); 235e334ad8fSJed Brown 236e334ad8fSJed Brown // Place PETSc vectors in CEED vectors 237e334ad8fSJed Brown ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr); 238e334ad8fSJed Brown ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr); 239e334ad8fSJed Brown CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, 240e334ad8fSJed Brown (PetscScalar *)q); 241e334ad8fSJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 242e334ad8fSJed Brown 243e334ad8fSJed Brown // Apply CEED operator 244e334ad8fSJed Brown CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, 245e334ad8fSJed Brown CEED_REQUEST_IMMEDIATE); 246e334ad8fSJed Brown 247e334ad8fSJed Brown // Restore vectors 248e334ad8fSJed Brown CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 249e334ad8fSJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 250e334ad8fSJed Brown ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr); 251e334ad8fSJed Brown ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr); 252e334ad8fSJed Brown 253e334ad8fSJed Brown // Local-to-Global 254e334ad8fSJed Brown ierr = VecZeroEntries(G); CHKERRQ(ierr); 255e334ad8fSJed Brown ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr); 256e334ad8fSJed Brown 257e334ad8fSJed Brown // Restore vectors 258e334ad8fSJed Brown ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr); 259e334ad8fSJed Brown PetscFunctionReturn(0); 260e334ad8fSJed Brown } 261e334ad8fSJed Brown 262e334ad8fSJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) { 263e334ad8fSJed Brown User user; 264e334ad8fSJed Brown Vec D_loc; 265e334ad8fSJed Brown PetscScalar *d; 266e334ad8fSJed Brown PetscMemType mem_type; 267e334ad8fSJed Brown 268e334ad8fSJed Brown PetscFunctionBeginUser; 269e334ad8fSJed Brown MatShellGetContext(A, &user); 270e334ad8fSJed Brown PetscCall(DMGetLocalVector(user->dm, &D_loc)); 271e334ad8fSJed Brown PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type)); 272e334ad8fSJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d); 273e334ad8fSJed Brown CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, 274e334ad8fSJed Brown CEED_REQUEST_IMMEDIATE); 275e334ad8fSJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL); 276e334ad8fSJed Brown PetscCall(VecRestoreArrayAndMemType(D_loc, &d)); 277e334ad8fSJed Brown PetscCall(VecZeroEntries(D)); 278e334ad8fSJed Brown PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D)); 279e334ad8fSJed Brown PetscCall(DMRestoreLocalVector(user->dm, &D_loc)); 280e334ad8fSJed Brown VecViewFromOptions(D, NULL, "-diag_vec_view"); 281e334ad8fSJed Brown PetscFunctionReturn(0); 282e334ad8fSJed Brown } 283e334ad8fSJed Brown 284544be873SJed Brown static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, 285544be873SJed Brown CeedVector *coo_values) { 286544be873SJed Brown PetscCount ncoo; 287544be873SJed Brown PetscInt *rows, *cols; 288544be873SJed Brown 289544be873SJed Brown PetscFunctionBeginUser; 290544be873SJed Brown if (pbdiagonal) { 291544be873SJed Brown CeedSize l_size; 292544be873SJed Brown CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL); 293544be873SJed Brown ncoo = l_size * 5; 294544be873SJed Brown rows = malloc(ncoo*sizeof(rows[0])); 295544be873SJed Brown cols = malloc(ncoo*sizeof(cols[0])); 296544be873SJed Brown for (PetscCount n=0; n<l_size/5; n++) { 297544be873SJed Brown for (PetscInt i=0; i<5; i++) { 298544be873SJed Brown for (PetscInt j=0; j<5; j++) { 299544be873SJed Brown rows[(n*5+i)*5+j] = n * 5 + i; 300544be873SJed Brown cols[(n*5+i)*5+j] = n * 5 + j; 301544be873SJed Brown } 302544be873SJed Brown } 303544be873SJed Brown } 304544be873SJed Brown } else { 305544be873SJed Brown PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, 306544be873SJed Brown &cols)); 307544be873SJed Brown } 308544be873SJed Brown PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols)); 309544be873SJed Brown free(rows); 310544be873SJed Brown free(cols); 311544be873SJed Brown CeedVectorCreate(user->ceed, ncoo, coo_values); 312544be873SJed Brown PetscFunctionReturn(0); 313544be873SJed Brown } 314544be873SJed Brown 315544be873SJed Brown static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, 316544be873SJed Brown CeedVector coo_values) { 317544be873SJed Brown CeedMemType mem_type = CEED_MEM_HOST; 318544be873SJed Brown const PetscScalar *values; 319544be873SJed Brown MatType mat_type; 320544be873SJed Brown 321544be873SJed Brown PetscFunctionBeginUser; 322544be873SJed Brown PetscCall(MatGetType(J, &mat_type)); 323544be873SJed Brown if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) 324544be873SJed Brown mem_type = CEED_MEM_DEVICE; 325544be873SJed Brown if (user->app_ctx->pmat_pbdiagonal) { 326544be873SJed Brown CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, 327544be873SJed Brown coo_values, CEED_REQUEST_IMMEDIATE); 328544be873SJed Brown } else { 329544be873SJed Brown CeedOperatorLinearAssemble(user->op_ijacobian, coo_values); 330544be873SJed Brown } 331544be873SJed Brown CeedVectorGetArrayRead(coo_values, mem_type, &values); 332544be873SJed Brown PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 333544be873SJed Brown CeedVectorRestoreArrayRead(coo_values, &values); 334544be873SJed Brown PetscFunctionReturn(0); 335544be873SJed Brown } 336544be873SJed Brown 337e334ad8fSJed Brown PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, 338e334ad8fSJed Brown PetscReal shift, Mat J, Mat J_pre, 339e334ad8fSJed Brown void *user_data) { 340e334ad8fSJed Brown User user = *(User *)user_data; 341d6bf345cSJed Brown PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 342e334ad8fSJed Brown PetscFunctionBeginUser; 343e334ad8fSJed Brown if (user->phys->ijacobian_time_shift_label) 344e334ad8fSJed Brown CeedOperatorContextSetDouble(user->op_ijacobian, 345e334ad8fSJed Brown user->phys->ijacobian_time_shift_label, &shift); 346d6bf345cSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 347e334ad8fSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 348e334ad8fSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, 349e334ad8fSJed Brown &J_pre_is_shell)); 350e334ad8fSJed Brown if (!user->matrices_set_up) { 351e334ad8fSJed Brown if (J_is_shell) { 352e334ad8fSJed Brown PetscCall(MatShellSetContext(J, user)); 353e334ad8fSJed Brown PetscCall(MatShellSetOperation(J, MATOP_MULT, 354e334ad8fSJed Brown (void (*)(void))MatMult_NS_IJacobian)); 355e334ad8fSJed Brown PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, 356e334ad8fSJed Brown (void (*)(void))MatGetDiagonal_NS_IJacobian)); 357e334ad8fSJed Brown PetscCall(MatSetUp(J)); 358e334ad8fSJed Brown } 359e334ad8fSJed Brown if (!J_pre_is_shell) { 360544be873SJed Brown PetscCall(FormPreallocation(user,user->app_ctx->pmat_pbdiagonal,J_pre, 361544be873SJed Brown &user->coo_values_pmat)); 362544be873SJed Brown } 363d6bf345cSJed Brown if (J != J_pre && !J_is_shell && !J_is_mffd) { 364544be873SJed Brown PetscCall(FormPreallocation(user,PETSC_FALSE,J, &user->coo_values_amat)); 365544be873SJed Brown } 366e334ad8fSJed Brown user->matrices_set_up = true; 367e334ad8fSJed Brown } 368e334ad8fSJed Brown if (!J_pre_is_shell) { 369544be873SJed Brown PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, 370544be873SJed Brown user->coo_values_pmat)); 371e334ad8fSJed Brown } 372d6bf345cSJed Brown if (user->coo_values_amat) { 373d6bf345cSJed Brown PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 374d6bf345cSJed Brown } else if (J_is_mffd) { 375d6bf345cSJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 376d6bf345cSJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 377d6bf345cSJed Brown } 378e334ad8fSJed Brown PetscFunctionReturn(0); 379e334ad8fSJed Brown } 380e334ad8fSJed Brown 381*f14660a4SJames Wright PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, 382*f14660a4SJames Wright PetscScalar time) { 38377841947SLeila Ghaffari Vec Q_loc; 38477841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 38577841947SLeila Ghaffari PetscViewer viewer; 38677841947SLeila Ghaffari PetscFunctionBeginUser; 38777841947SLeila Ghaffari 38877841947SLeila Ghaffari // Set up output 389*f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 390*f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 391*f14660a4SJames Wright PetscCall(VecZeroEntries(Q_loc)); 392*f14660a4SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 39377841947SLeila Ghaffari 39477841947SLeila Ghaffari // Output 395*f14660a4SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, 39608140895SJed Brown "%s/ns-%03" PetscInt_FMT ".vtu", 397*f14660a4SJames Wright user->app_ctx->output_dir, step_no + user->app_ctx->cont_steps)); 398*f14660a4SJames Wright 399*f14660a4SJames Wright PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, 400*f14660a4SJames Wright FILE_MODE_WRITE, &viewer)); 401*f14660a4SJames Wright PetscCall(VecView(Q_loc, viewer)); 402*f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 40377841947SLeila Ghaffari if (user->dm_viz) { 40477841947SLeila Ghaffari Vec Q_refined, Q_refined_loc; 40577841947SLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 40677841947SLeila Ghaffari PetscViewer viewer_refined; 40777841947SLeila Ghaffari 408*f14660a4SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 409*f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 410*f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 411*f14660a4SJames Wright 412*f14660a4SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 413*f14660a4SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 414*f14660a4SJames Wright PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, 415*f14660a4SJames Wright Q_refined_loc)); 416*f14660a4SJames Wright 417*f14660a4SJames Wright PetscCall(PetscSNPrintf(file_path_refined, sizeof file_path_refined, 41808140895SJed Brown "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, 419*f14660a4SJames Wright step_no + user->app_ctx->cont_steps)); 420*f14660a4SJames Wright 421*f14660a4SJames Wright PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), 422*f14660a4SJames Wright file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 423*f14660a4SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 424*f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 425*f14660a4SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 426*f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 42777841947SLeila Ghaffari } 428*f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 42977841947SLeila Ghaffari 43077841947SLeila Ghaffari // Save data in a binary file for continuation of simulations 431*f14660a4SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", 432*f14660a4SJames Wright user->app_ctx->output_dir)); 433*f14660a4SJames Wright PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, 434*f14660a4SJames Wright &viewer)); 435*f14660a4SJames Wright 436*f14660a4SJames Wright PetscCall(VecView(Q, viewer)); 437*f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 43877841947SLeila Ghaffari 43977841947SLeila Ghaffari // Save time stamp 44077841947SLeila Ghaffari // Dimensionalize time back 44177841947SLeila Ghaffari time /= user->units->second; 442*f14660a4SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin", 443*f14660a4SJames Wright user->app_ctx->output_dir)); 444*f14660a4SJames Wright PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, 445*f14660a4SJames Wright &viewer)); 446*f14660a4SJames Wright 44777841947SLeila Ghaffari #if PETSC_VERSION_GE(3,13,0) 448*f14660a4SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 44977841947SLeila Ghaffari #else 450*f14660a4SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true)); 45177841947SLeila Ghaffari #endif 452*f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 453*f14660a4SJames Wright 454*f14660a4SJames Wright PetscFunctionReturn(0); 455*f14660a4SJames Wright } 456*f14660a4SJames Wright 457*f14660a4SJames Wright // User provided TS Monitor 458*f14660a4SJames Wright PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, 459*f14660a4SJames Wright Vec Q, void *ctx) { 460*f14660a4SJames Wright User user = ctx; 461*f14660a4SJames Wright PetscFunctionBeginUser; 462*f14660a4SJames Wright 463*f14660a4SJames Wright // Print every 'output_freq' steps 464*f14660a4SJames Wright if (user->app_ctx->output_freq <= 0 465*f14660a4SJames Wright || step_no % user->app_ctx->output_freq != 0) 466*f14660a4SJames Wright PetscFunctionReturn(0); 467*f14660a4SJames Wright 468*f14660a4SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 46977841947SLeila Ghaffari 47077841947SLeila Ghaffari PetscFunctionReturn(0); 47177841947SLeila Ghaffari } 47277841947SLeila Ghaffari 47377841947SLeila Ghaffari // TS: Create, setup, and solve 47477841947SLeila Ghaffari PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, 47577841947SLeila Ghaffari Vec *Q, PetscScalar *f_time, TS *ts) { 47677841947SLeila Ghaffari MPI_Comm comm = user->comm; 47777841947SLeila Ghaffari TSAdapt adapt; 47877841947SLeila Ghaffari PetscScalar final_time; 47977841947SLeila Ghaffari PetscErrorCode ierr; 48077841947SLeila Ghaffari PetscFunctionBeginUser; 48177841947SLeila Ghaffari 48277841947SLeila Ghaffari ierr = TSCreate(comm, ts); CHKERRQ(ierr); 48377841947SLeila Ghaffari ierr = TSSetDM(*ts, dm); CHKERRQ(ierr); 48477841947SLeila Ghaffari if (phys->implicit) { 48577841947SLeila Ghaffari ierr = TSSetType(*ts, TSBDF); CHKERRQ(ierr); 48677841947SLeila Ghaffari if (user->op_ifunction) { 48777841947SLeila Ghaffari ierr = TSSetIFunction(*ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 48877841947SLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 48977841947SLeila Ghaffari ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 49077841947SLeila Ghaffari } 491e334ad8fSJed Brown if (user->op_ijacobian) { 492e334ad8fSJed Brown ierr = DMTSSetIJacobian(dm, FormIJacobian_NS, &user); CHKERRQ(ierr); 493544be873SJed Brown if (app_ctx->amat_type) { 494544be873SJed Brown Mat Pmat,Amat; 495544be873SJed Brown ierr = DMCreateMatrix(dm, &Pmat); CHKERRQ(ierr); 496544be873SJed Brown ierr = DMSetMatType(dm, app_ctx->amat_type); CHKERRQ(ierr); 497544be873SJed Brown ierr = DMCreateMatrix(dm, &Amat); CHKERRQ(ierr); 498544be873SJed Brown ierr = TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL); CHKERRQ(ierr); 499544be873SJed Brown ierr = MatDestroy(&Amat); CHKERRQ(ierr); 500544be873SJed Brown ierr = MatDestroy(&Pmat); CHKERRQ(ierr); 501544be873SJed Brown } 502e334ad8fSJed Brown } 50377841947SLeila Ghaffari } else { 50477841947SLeila Ghaffari if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, 50577841947SLeila Ghaffari "Problem does not provide RHSFunction"); 50677841947SLeila Ghaffari ierr = TSSetType(*ts, TSRK); CHKERRQ(ierr); 50777841947SLeila Ghaffari ierr = TSRKSetType(*ts, TSRK5F); CHKERRQ(ierr); 50877841947SLeila Ghaffari ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 50977841947SLeila Ghaffari } 51077841947SLeila Ghaffari ierr = TSSetMaxTime(*ts, 500. * user->units->second); CHKERRQ(ierr); 51177841947SLeila Ghaffari ierr = TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 51277841947SLeila Ghaffari ierr = TSSetTimeStep(*ts, 1.e-2 * user->units->second); CHKERRQ(ierr); 51377841947SLeila Ghaffari if (app_ctx->test_mode) {ierr = TSSetMaxSteps(*ts, 10); CHKERRQ(ierr);} 51477841947SLeila Ghaffari ierr = TSGetAdapt(*ts, &adapt); CHKERRQ(ierr); 51577841947SLeila Ghaffari ierr = TSAdaptSetStepLimits(adapt, 51677841947SLeila Ghaffari 1.e-12 * user->units->second, 51777841947SLeila Ghaffari 1.e2 * user->units->second); CHKERRQ(ierr); 51877841947SLeila Ghaffari ierr = TSSetFromOptions(*ts); CHKERRQ(ierr); 5195e82a6e1SJeremy L Thompson user->time = -1.0; // require all BCs and ctx to be updated 5205e82a6e1SJeremy L Thompson user->dt = -1.0; 52177841947SLeila Ghaffari if (!app_ctx->cont_steps) { // print initial condition 52277841947SLeila Ghaffari if (!app_ctx->test_mode) { 52377841947SLeila Ghaffari ierr = TSMonitor_NS(*ts, 0, 0., *Q, user); CHKERRQ(ierr); 52477841947SLeila Ghaffari } 52577841947SLeila Ghaffari } else { // continue from time of last output 52677841947SLeila Ghaffari PetscReal time; 52777841947SLeila Ghaffari PetscInt count; 52877841947SLeila Ghaffari PetscViewer viewer; 52977841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 53077841947SLeila Ghaffari ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin", 53177841947SLeila Ghaffari app_ctx->output_dir); CHKERRQ(ierr); 53277841947SLeila Ghaffari ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer); 53377841947SLeila Ghaffari CHKERRQ(ierr); 53477841947SLeila Ghaffari ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 53577841947SLeila Ghaffari CHKERRQ(ierr); 53677841947SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 53777841947SLeila Ghaffari ierr = TSSetTime(*ts, time * user->units->second); CHKERRQ(ierr); 53877841947SLeila Ghaffari } 53977841947SLeila Ghaffari if (!app_ctx->test_mode) { 54077841947SLeila Ghaffari ierr = TSMonitorSet(*ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 54177841947SLeila Ghaffari } 54277841947SLeila Ghaffari 54377841947SLeila Ghaffari // Solve 5447b39487dSJeremy L Thompson PetscScalar start_time; 5457b39487dSJeremy L Thompson ierr = TSGetTime(*ts, &start_time); CHKERRQ(ierr); 5467b39487dSJeremy L Thompson 5477b39487dSJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 5487b39487dSJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 5497b39487dSJeremy L Thompson PetscCall(TSSetStepNumber(*ts, 0)); 5507b39487dSJeremy L Thompson if (PetscPreLoadingOn) { 5517b39487dSJeremy L Thompson // LCOV_EXCL_START 5527b39487dSJeremy L Thompson SNES snes; 5537b39487dSJeremy L Thompson Vec Q_preload; 5547b39487dSJeremy L Thompson PetscReal rtol; 5557b39487dSJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 5567b39487dSJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 5577b39487dSJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 5587b39487dSJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 5597b39487dSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, 5607b39487dSJeremy L Thompson PETSC_DEFAULT, PETSC_DEFAULT)); 5617b39487dSJeremy L Thompson PetscCall(TSSetSolution(*ts, *Q)); 5627b39487dSJeremy L Thompson PetscCall(TSStep(*ts)); 5637b39487dSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, 5647b39487dSJeremy L Thompson PETSC_DEFAULT, PETSC_DEFAULT)); 5657b39487dSJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 5667b39487dSJeremy L Thompson // LCOV_EXCL_STOP 5677b39487dSJeremy L Thompson } else { 56877841947SLeila Ghaffari ierr = PetscBarrier((PetscObject) *ts); CHKERRQ(ierr); 56977841947SLeila Ghaffari ierr = TSSolve(*ts, *Q); CHKERRQ(ierr); 5707b39487dSJeremy L Thompson } 5717b39487dSJeremy L Thompson PetscPreLoadEnd(); 5727b39487dSJeremy L Thompson 5737b39487dSJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 57477841947SLeila Ghaffari *f_time = final_time; 5757b39487dSJeremy L Thompson 576*f14660a4SJames Wright 57777841947SLeila Ghaffari if (!app_ctx->test_mode) { 578*f14660a4SJames Wright if (user->app_ctx->output_freq > 0 || user->app_ctx->output_freq == -1) { 579*f14660a4SJames Wright PetscInt step_no; 580*f14660a4SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 581*f14660a4SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 582*f14660a4SJames Wright } 583*f14660a4SJames Wright 5847b39487dSJeremy L Thompson PetscLogEvent stage_id; 5857b39487dSJeremy L Thompson PetscStageLog stage_log; 5867b39487dSJeremy L Thompson 5877b39487dSJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 5887b39487dSJeremy L Thompson PetscCall(PetscLogGetStageLog(&stage_log)); 5897b39487dSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, 59077841947SLeila Ghaffari "Time taken for solution (sec): %g\n", 5917b39487dSJeremy L Thompson stage_log->stageInfo[stage_id].perfInfo.time)); 59277841947SLeila Ghaffari } 59377841947SLeila Ghaffari PetscFunctionReturn(0); 59477841947SLeila Ghaffari } 595