xref: /libCEED/examples/fluids/src/setupts.c (revision 04292f7dad43a72a08e1eaa490df719405885952)
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 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
1349aac155SJeremy L Thompson #include <petscts.h>
1449aac155SJeremy L Thompson 
1577841947SLeila Ghaffari #include "../navierstokes.h"
16ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
1777841947SLeila Ghaffari 
18*04292f7dSJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
19*04292f7dSJames Wright PetscErrorCode CreateKspMassOperator(User user, CeedData ceed_data) {
20*04292f7dSJames Wright   Ceed                 ceed = user->ceed;
21*04292f7dSJames Wright   DM                   dm   = user->dm;
2277841947SLeila Ghaffari   CeedQFunction        qf_mass;
2377841947SLeila Ghaffari   CeedOperator         op_mass;
24*04292f7dSJames Wright   OperatorApplyContext mass_matop_ctx;
2577841947SLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
2677841947SLeila Ghaffari 
27f17d818dSJames Wright   PetscFunctionBeginUser;
28a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
29a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
3077841947SLeila Ghaffari 
31ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
32a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
33a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
34356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
35a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3677841947SLeila Ghaffari 
37*04292f7dSJames Wright   {  // -- Setup KSP for mass operator
38*04292f7dSJames Wright     Mat      mat_mass;
39*04292f7dSJames Wright     Vec      Ones_loc;
40*04292f7dSJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
4177841947SLeila Ghaffari 
42*04292f7dSJames Wright     PetscCall(DMCreateLocalVector(dm, &Ones_loc));
433221f4d3SJames Wright     PetscCall(VecSet(Ones_loc, 1));
44*04292f7dSJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_mass, NULL, NULL, Ones_loc, NULL, &mass_matop_ctx));
45*04292f7dSJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
4677841947SLeila Ghaffari 
47*04292f7dSJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
48*04292f7dSJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
49*04292f7dSJames Wright     {  // lumped by default
50*04292f7dSJames Wright       PC pc;
51*04292f7dSJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
52*04292f7dSJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
53*04292f7dSJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
54*04292f7dSJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
55*04292f7dSJames Wright     }
56*04292f7dSJames Wright     PetscCall(KSPSetOperators(user->mass_ksp, mat_mass, mat_mass));
57*04292f7dSJames Wright     PetscCall(KSPSetFromOptions(user->mass_ksp));
58*04292f7dSJames Wright     PetscCall(VecDestroy(&Ones_loc));
59*04292f7dSJames Wright   }
6077841947SLeila Ghaffari 
6177841947SLeila Ghaffari   // Cleanup
62a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
63a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
64ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6577841947SLeila Ghaffari }
6677841947SLeila Ghaffari 
67a0b9a424SJames Wright // Insert Boundary values if it's a new time
68a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
69a0b9a424SJames Wright   PetscFunctionBeginUser;
70a0b9a424SJames Wright   if (user->time_bc_set != t) {
71a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
72a0b9a424SJames Wright     user->time_bc_set = t;
73a0b9a424SJames Wright   }
74ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
75a0b9a424SJames Wright }
76a0b9a424SJames Wright 
7777841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
7877841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
7977841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
8077841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
8177841947SLeila Ghaffari   User        user = *(User *)user_data;
827d4c6defSJames Wright   Ceed        ceed = user->ceed;
83c798d55aSJames Wright   PetscScalar dt;
849ad5e8e4SJames Wright   Vec         Q_loc = user->Q_loc;
8577841947SLeila Ghaffari 
86f17d818dSJames Wright   PetscFunctionBeginUser;
875e82a6e1SJeremy L Thompson   // Update time dependent data
88a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
897d4c6defSJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
902b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
917d4c6defSJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
9277841947SLeila Ghaffari 
939ad5e8e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
9477841947SLeila Ghaffari 
95d6e67e47SJames Wright   // Inverse of the lumped mass matrix
96*04292f7dSJames Wright   PetscCall(KSPSolve(user->mass_ksp, G, G));
97ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
9877841947SLeila Ghaffari }
9977841947SLeila Ghaffari 
100ca69d878SAdeleke O. Bankole // Surface forces function setup
101ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
102ca69d878SAdeleke O. Bankole   DMLabel            face_label;
103ca69d878SAdeleke O. Bankole   const PetscScalar *g;
104d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
105ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
106d6734f85SAdeleke O. Bankole   PetscSection       s;
107ca69d878SAdeleke O. Bankole 
108ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
109ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
110ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
111ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
112ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
113ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
114ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
115ca69d878SAdeleke O. Bankole     IS             wall_is;
116d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
117ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
118ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
119ca69d878SAdeleke O. Bankole       PetscInt        num_points;
120d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
121ca69d878SAdeleke O. Bankole       const PetscInt *points;
122d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
123ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
124ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
125ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
126ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
127ca69d878SAdeleke O. Bankole         const StateConservative *r;
128ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
129d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
130d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
131ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
132d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
133d6734f85SAdeleke O. Bankole           }
134ca69d878SAdeleke O. Bankole         }
135ca69d878SAdeleke O. Bankole       }
136ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
137ca69d878SAdeleke O. Bankole     }
138ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
139ca69d878SAdeleke O. Bankole   }
140ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
141ca69d878SAdeleke O. Bankole   //  Restore Vectors
142ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
143ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
144ca69d878SAdeleke O. Bankole }
145ca69d878SAdeleke O. Bankole 
14677841947SLeila Ghaffari // Implicit time-stepper function setup
1472b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
14877841947SLeila Ghaffari   User         user = *(User *)user_data;
1497d4c6defSJames Wright   Ceed         ceed = user->ceed;
150c798d55aSJames Wright   PetscScalar  dt;
1515e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
15277841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
15377841947SLeila Ghaffari 
154f17d818dSJames Wright   PetscFunctionBeginUser;
1555e82a6e1SJeremy L Thompson   // Get local vectors
156ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1575e82a6e1SJeremy L Thompson 
1585e82a6e1SJeremy L Thompson   // Update time dependent data
159a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
1607d4c6defSJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1612b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
1627d4c6defSJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
16377841947SLeila Ghaffari 
16477841947SLeila Ghaffari   // Global-to-local
165878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
166878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
167878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
168878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
16977841947SLeila Ghaffari 
17077841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
171c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
172c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
173c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
17477841947SLeila Ghaffari 
17577841947SLeila Ghaffari   // Apply CEED operator
17675d1798cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
17775d1798cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
178a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
17975d1798cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
18075d1798cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
18177841947SLeila Ghaffari 
18277841947SLeila Ghaffari   // Restore vectors
183c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
184c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
185c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
18677841947SLeila Ghaffari 
18719ffbc25SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
188cbef7084SJames Wright     PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc));
18919ffbc25SJames Wright   }
190be34b3b7SJames Wright 
19177841947SLeila Ghaffari   // Local-to-Global
1922b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1932b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
19477841947SLeila Ghaffari 
19577841947SLeila Ghaffari   // Restore vectors
196ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
197ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19877841947SLeila Ghaffari }
19977841947SLeila Ghaffari 
2002b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
201544be873SJed Brown   PetscCount ncoo;
202457e73b2SJames Wright   PetscInt  *rows_petsc, *cols_petsc;
20301f0e615SJames Wright   CeedInt   *rows_ceed, *cols_ceed;
204544be873SJed Brown 
205544be873SJed Brown   PetscFunctionBeginUser;
206544be873SJed Brown   if (pbdiagonal) {
20701f0e615SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
208544be873SJed Brown   } else {
209a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
21001f0e615SJames Wright   }
211457e73b2SJames Wright   PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
212457e73b2SJames Wright   PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
213457e73b2SJames Wright   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
214457e73b2SJames Wright   free(rows_petsc);
215457e73b2SJames Wright   free(cols_petsc);
216a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values));
217ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
218544be873SJed Brown }
219544be873SJed Brown 
2202b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
221544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
222544be873SJed Brown   const PetscScalar *values;
223544be873SJed Brown   MatType            mat_type;
224544be873SJed Brown 
225544be873SJed Brown   PetscFunctionBeginUser;
226544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
2272b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
228fbb546ddSJames Wright   if (pbdiagonal) {
22975d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
23075d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
231a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE));
23275d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
23375d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
234544be873SJed Brown   } else {
23575d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
23675d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
237a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values));
23875d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
23975d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
240544be873SJed Brown   }
241a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values));
242544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
243a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values));
244ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
245544be873SJed Brown }
246544be873SJed Brown 
2472b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
248e334ad8fSJed Brown   User      user = *(User *)user_data;
249a424bcd0SJames Wright   Ceed      ceed = user->ceed;
250d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
251f17d818dSJames Wright 
252e334ad8fSJed Brown   PetscFunctionBeginUser;
253a424bcd0SJames Wright   if (user->phys->ijacobian_time_shift_label)
254a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
255d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
256e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2572b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
258e334ad8fSJed Brown   if (!user->matrices_set_up) {
259e334ad8fSJed Brown     if (J_is_shell) {
260b07979f9SJames Wright       OperatorApplyContext op_ijacobian_ctx;
261b07979f9SJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
262b07979f9SJames Wright                                  &op_ijacobian_ctx);
263b07979f9SJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
264b07979f9SJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
265b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
266b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
267e334ad8fSJed Brown       PetscCall(MatSetUp(J));
268e334ad8fSJed Brown     }
269e334ad8fSJed Brown     if (!J_pre_is_shell) {
2702b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
271544be873SJed Brown     }
272d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
273544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
274544be873SJed Brown     }
275e334ad8fSJed Brown     user->matrices_set_up = true;
276e334ad8fSJed Brown   }
277e334ad8fSJed Brown   if (!J_pre_is_shell) {
2782b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
279e334ad8fSJed Brown   }
280d6bf345cSJed Brown   if (user->coo_values_amat) {
281d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
282d6bf345cSJed Brown   } else if (J_is_mffd) {
283d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
284d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
285d6bf345cSJed Brown   }
286ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
287e334ad8fSJed Brown }
288e334ad8fSJed Brown 
2892b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
29077841947SLeila Ghaffari   Vec         Q_loc;
29177841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
29277841947SLeila Ghaffari   PetscViewer viewer;
29377841947SLeila Ghaffari 
294f17d818dSJames Wright   PetscFunctionBeginUser;
29537cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
29677841947SLeila Ghaffari     // Set up output
297f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
298f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
299f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
300f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
30177841947SLeila Ghaffari 
30277841947SLeila Ghaffari     // Output
30337cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
304f14660a4SJames Wright 
3052b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
306f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
307f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
30877841947SLeila Ghaffari     if (user->dm_viz) {
30977841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
31077841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
31177841947SLeila Ghaffari       PetscViewer viewer_refined;
31277841947SLeila Ghaffari 
313f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
314f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
315f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
316f14660a4SJames Wright 
317f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
318f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3192b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
320f14660a4SJames Wright 
32137cbb16aSJed Brown       PetscCall(
32237cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
323f14660a4SJames Wright 
3242b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
325f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
326f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
327f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
328f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
32977841947SLeila Ghaffari     }
330f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
33137cbb16aSJed Brown   }
33277841947SLeila Ghaffari 
33377841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
33469293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
33537cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
33669293791SJames Wright   } else {
3372b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
33869293791SJames Wright   }
3392b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
340f14660a4SJames Wright 
3410de6a49fSJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
3420de6a49fSJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
3434de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3444de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
3454de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
346f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
347f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
348ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
349f14660a4SJames Wright }
350f14660a4SJames Wright 
351ca69d878SAdeleke O. Bankole // CSV Monitor
352ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
353ca69d878SAdeleke O. Bankole   User              user = ctx;
354ca69d878SAdeleke O. Bankole   Vec               G_loc;
355ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
356ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
357ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
358ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
359ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
360ca69d878SAdeleke O. Bankole   PetscBool         iascii;
361ca69d878SAdeleke O. Bankole 
362ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
363ee4ca9cbSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
364ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
365ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
366ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
367ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
368ca69d878SAdeleke O. Bankole 
369ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
370ca69d878SAdeleke O. Bankole 
371ca69d878SAdeleke O. Bankole   if (iascii) {
372ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
373ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
374ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
375ca69d878SAdeleke O. Bankole     }
376ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
377ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
378ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
379ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
380ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
381ca69d878SAdeleke O. Bankole 
382ca69d878SAdeleke O. Bankole       } else {
383ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
384ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
385ca69d878SAdeleke O. Bankole       }
386ca69d878SAdeleke O. Bankole     }
387ca69d878SAdeleke O. Bankole   }
388ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
389ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
390ca69d878SAdeleke O. Bankole }
391ca69d878SAdeleke O. Bankole 
392f14660a4SJames Wright // User provided TS Monitor
3932b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
394f14660a4SJames Wright   User user = ctx;
395f14660a4SJames Wright 
396f17d818dSJames Wright   PetscFunctionBeginUser;
39737cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
398894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
39949aac155SJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
400ee4ca9cbSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
40149aac155SJeremy L Thompson   }
402f14660a4SJames Wright 
403f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
404ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
40577841947SLeila Ghaffari }
40677841947SLeila Ghaffari 
40777841947SLeila Ghaffari // TS: Create, setup, and solve
4082b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
40977841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
41077841947SLeila Ghaffari   TSAdapt     adapt;
41177841947SLeila Ghaffari   PetscScalar final_time;
41277841947SLeila Ghaffari 
413f17d818dSJames Wright   PetscFunctionBeginUser;
4142b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4152b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
4161a7db67cSJames Wright   PetscCall(TSSetApplicationContext(*ts, user));
41777841947SLeila Ghaffari   if (phys->implicit) {
4182b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
41977841947SLeila Ghaffari     if (user->op_ifunction) {
4202b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
42177841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4222b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
42377841947SLeila Ghaffari     }
424e334ad8fSJed Brown     if (user->op_ijacobian) {
4252b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
426544be873SJed Brown       if (app_ctx->amat_type) {
427544be873SJed Brown         Mat Pmat, Amat;
4282b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4292b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4302b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4312b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4322b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4332b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
434544be873SJed Brown       }
435e334ad8fSJed Brown     }
43677841947SLeila Ghaffari   } else {
4379ad5e8e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4382b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4392b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4402b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
44177841947SLeila Ghaffari   }
4422b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4432b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
44451b00d91SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4452b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4462b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4472b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4482b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
44915637395SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
450d55646a4SJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
45177841947SLeila Ghaffari     PetscInt    count;
45277841947SLeila Ghaffari     PetscViewer viewer;
4532b730f8bSJeremy L Thompson 
4544de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4552b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4564de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4572b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4584de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4594de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4604de8550aSJed Brown     }
4614de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
462d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
46377841947SLeila Ghaffari   }
4648fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4652b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4668fb33541SJames Wright   }
467ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
468ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
469ca69d878SAdeleke O. Bankole   }
470b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
471f5452247SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
472495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
473a424bcd0SJames Wright     PetscCallCeed(user->ceed,
474a424bcd0SJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
475a175e481SJames Wright   }
4764e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
47777841947SLeila Ghaffari 
47828134cfaSJames Wright   if (app_ctx->sgs_train_enable) {
47928134cfaSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
48028134cfaSJames Wright     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
48128134cfaSJames Wright   }
48277841947SLeila Ghaffari   // Solve
483d8e0aecdSJed Brown   PetscReal start_time;
484d8e0aecdSJed Brown   PetscInt  start_step;
4852b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
486d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
4877b39487dSJeremy L Thompson 
488711a423aSJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
4897b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
4907b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
491d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
4927b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
4937b39487dSJeremy L Thompson     // LCOV_EXCL_START
4947b39487dSJeremy L Thompson     SNES      snes;
4957b39487dSJeremy L Thompson     Vec       Q_preload;
4967b39487dSJeremy L Thompson     PetscReal rtol;
4977b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
4987b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
4997b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5007b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5012b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
502fc4c3d69SJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
5037b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
5042b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5057b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
5067b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
5077b39487dSJeremy L Thompson   } else {
5082b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5092b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
5107b39487dSJeremy L Thompson   }
5117b39487dSJeremy L Thompson   PetscPreLoadEnd();
5127b39487dSJeremy L Thompson 
5137b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
51477841947SLeila Ghaffari   *f_time = final_time;
5157b39487dSJeremy L Thompson 
5168fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
517f14660a4SJames Wright     PetscInt step_no;
518f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
519a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
520f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
521f14660a4SJames Wright     }
522f14660a4SJames Wright 
52375d1798cSJames Wright     PetscLogStage      stage_id;
524711a423aSJed Brown     PetscEventPerfInfo stage_perf;
5257b39487dSJeremy L Thompson 
5267b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
527711a423aSJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
528711a423aSJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
52977841947SLeila Ghaffari   }
530ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53177841947SLeila Ghaffari }
532