1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 #include "../qfunctions/monitor_totalkineticenergy.h" 5 #include <navierstokes.h> 6 #include <time.h> 7 8 typedef struct { 9 OperatorApplyContext op_monitor_ctx; 10 Vec values; 11 PetscScalar *sum_values; 12 PetscScalar volume; 13 PetscInt num_comps; 14 PetscBool is_header_written; 15 } *MonitorDissipation; 16 17 PetscErrorCode MonitorDissipationDestroy(void **ctx) { 18 MonitorDissipation monitor_ctx = *(MonitorDissipation *)ctx; 19 20 PetscFunctionBeginUser; 21 PetscCall(OperatorApplyContextDestroy(monitor_ctx->op_monitor_ctx)); 22 PetscCall(VecDestroy(&monitor_ctx->values)); 23 PetscCall(PetscFree(monitor_ctx->sum_values)); 24 PetscCall(PetscFree(monitor_ctx)); 25 PetscFunctionReturn(PETSC_SUCCESS); 26 } 27 28 /** @brief Create CeedElemRestriction for collocated data in component-major order. 29 a. Sets the strides of the restriction to component-major order 30 Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction. 31 */ 32 static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 33 CeedElemRestriction *elem_restr_collocated) { 34 CeedInt num_elem_qpts, loc_num_elem; 35 36 PetscFunctionBeginUser; 37 PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts)); 38 PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem)); 39 40 const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 41 PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 42 elem_restr_collocated)); 43 PetscFunctionReturn(PETSC_SUCCESS); 44 } 45 46 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx) { 47 Honee honee; 48 Ceed ceed; 49 CeedQFunction qf_monitor; 50 CeedOperator op_monitor; 51 CeedElemRestriction elem_restr_qd, elem_restr_totalke, elem_restr_q; 52 CeedBasis basis_q; 53 CeedVector q_data; 54 CeedInt num_comp_q, q_data_size; 55 DMLabel domain_label = NULL; 56 PetscInt label_value = 0, num_comp_totalke = 5, dim; 57 MonitorDissipation monitor_ctx; 58 CeedQFunctionContext newt_qfctx; 59 PetscBool is_ascii; 60 61 PetscFunctionBeginUser; 62 PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii)); 63 PetscCheck(is_ascii, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only supports ASCII viewers"); 64 PetscCall(TSGetApplicationContext(ts, &honee)); 65 PetscCall(DMGetDimension(honee->dm, &dim)); 66 ceed = honee->ceed; 67 68 PetscCall(PetscNew(&monitor_ctx)); 69 PetscCall(HoneeCalculateDomainSize(honee, &monitor_ctx->volume)); 70 71 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 72 PetscCall(QDataGet(ceed, honee->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 73 &q_data_size)); 74 75 { // Get restriction and basis from the RHS function 76 CeedOperator *sub_ops; 77 CeedOperatorField op_field; 78 PetscInt sub_op_index = 0; // will be 0 for the volume op 79 80 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 81 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 82 83 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 84 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newt_qfctx)); 85 } 86 PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_totalke, basis_q, elem_restr_q, &elem_restr_totalke)); 87 88 switch (honee->phys->state_var) { 89 case STATEVAR_PRIMITIVE: 90 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Prim, MonitorTotalKineticEnergy_Prim_loc, &qf_monitor)); 91 break; 92 case STATEVAR_CONSERVATIVE: 93 PetscCallCeed(ceed, 94 CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Conserv, MonitorTotalKineticEnergy_Conserv_loc, &qf_monitor)); 95 break; 96 case STATEVAR_ENTROPY: 97 PetscCallCeed(ceed, 98 CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Entropy, MonitorTotalKineticEnergy_Entropy_loc, &qf_monitor)); 99 break; 100 } 101 102 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_monitor, newt_qfctx)); 103 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q", num_comp_q, CEED_EVAL_INTERP)); 104 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 105 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q_data", q_data_size, CEED_EVAL_NONE)); 106 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_monitor, "v", num_comp_totalke, CEED_EVAL_NONE)); 107 108 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_monitor, NULL, NULL, &op_monitor)); 109 PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 110 PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 111 PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data)); 112 PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "v", elem_restr_totalke, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 113 114 PetscCall(OperatorApplyContextCreate(honee->dm, NULL, ceed, op_monitor, NULL, NULL, NULL, NULL, &monitor_ctx->op_monitor_ctx)); 115 116 PetscCall(CeedOperatorCreateLocalVecs(op_monitor, DMReturnVecType(honee->dm), PETSC_COMM_SELF, NULL, &monitor_ctx->values)); 117 PetscCall(VecSetBlockSize(monitor_ctx->values, num_comp_totalke)); 118 monitor_ctx->num_comps = num_comp_totalke; 119 PetscCall(PetscMalloc1(num_comp_totalke, &monitor_ctx->sum_values)); 120 121 ctx->data = monitor_ctx; 122 ctx->data_destroy = MonitorDissipationDestroy; 123 124 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newt_qfctx)); 125 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 126 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 127 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_totalke)); 128 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 129 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 130 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_monitor)); 131 PetscCallCeed(ceed, CeedOperatorDestroy(&op_monitor)); 132 PetscFunctionReturn(PETSC_SUCCESS); 133 } 134 135 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) { 136 MonitorDissipation monitor_ctx = (MonitorDissipation)ctx->data; 137 Honee honee; 138 MPI_Comm comm; 139 PetscMPIInt rank; 140 TSConvergedReason reason; 141 static const char *field_names[] = {"TotalKineticEnergy", "StrainDissipation", "DivergenceDissipation", "VolumeExpansion", "MuVorticity2"}; 142 143 PetscFunctionBeginUser; 144 PetscCall(TSGetConvergedReason(ts, &reason)); 145 if (!(step % ctx->view_interval == 0 || reason != TS_CONVERGED_ITERATING)) PetscFunctionReturn(PETSC_SUCCESS); 146 PetscCall(TSGetApplicationContext(ts, &honee)); 147 PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 148 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 149 150 PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time)); 151 PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc)); 152 153 PetscCall(ApplyCeedOperatorLocalToLocal(honee->Q_loc, monitor_ctx->values, monitor_ctx->op_monitor_ctx)); 154 155 for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) { 156 PetscCall(VecStrideSum(monitor_ctx->values, i, &monitor_ctx->sum_values[i])); 157 } 158 159 if (rank == 0) PetscCallMPI(MPI_Reduce(MPI_IN_PLACE, monitor_ctx->sum_values, monitor_ctx->num_comps, MPIU_SCALAR, MPI_SUM, 0, comm)); 160 else PetscCallMPI(MPI_Reduce(monitor_ctx->sum_values, monitor_ctx->sum_values, monitor_ctx->num_comps, MPIU_SCALAR, MPI_SUM, 0, comm)); 161 162 if (rank == 0) 163 for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) monitor_ctx->sum_values[i] /= monitor_ctx->volume; 164 165 if (ctx->format == PETSC_VIEWER_ASCII_CSV) { 166 if (!monitor_ctx->is_header_written) { 167 char buf[PETSC_MAX_PATH_LEN]; 168 const char *buf_const; 169 time_t t = time(NULL); 170 171 PetscCall(PetscGetVersion(buf, sizeof(buf))); 172 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_version: %s\n", buf)); 173 PetscCall(PetscGetPetscDir(&buf_const)); 174 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_directory: %s\n", buf_const)); 175 PetscCall(PetscGetArchType(buf, sizeof(buf))); 176 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_arch: %s\n", buf)); 177 if (strftime(buf, sizeof(buf), "%FT%T%z", localtime(&t)) == 0) SETERRQ(comm, PETSC_ERR_SYS, "strftime() call failed"); 178 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date: %s\n", buf)); 179 if (strftime(buf, sizeof(buf), "%Z", localtime(&t)) == 0) SETERRQ(comm, PETSC_ERR_SYS, "strftime() call failed"); 180 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date_timezone: %s\n", buf)); 181 PetscCall(PetscGetUserName(buf, sizeof(buf))); 182 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# username: %s\n", buf)); 183 PetscCall(PetscGetHostName(buf, sizeof(buf))); 184 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# hostname: %s\n", buf)); 185 PetscCall(PetscGetWorkingDirectory(buf, sizeof(buf))); 186 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# working_directory: %s\n", buf)); 187 188 PetscCall(PetscViewerFileGetName(ctx->viewer, &buf_const)); 189 PetscCall(PetscGetFullPath(buf_const, buf, sizeof(buf))); 190 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# original_file_path: %s\n", buf)); 191 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "#\n")); 192 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "step,time,")); 193 for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) { 194 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%s", field_names[i])); 195 if (i < monitor_ctx->num_comps - 1) PetscCall(PetscViewerASCIIPrintf(ctx->viewer, ",")); 196 } 197 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n")); 198 monitor_ctx->is_header_written = PETSC_TRUE; 199 } 200 201 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%" PetscInt_FMT ",%0.17e,", step, solution_time)); 202 for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) { 203 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%0.17e", monitor_ctx->sum_values[i])); 204 if (i < monitor_ctx->num_comps - 1) PetscCall(PetscViewerASCIIPrintf(ctx->viewer, ",")); 205 } 206 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n")); 207 } else { 208 for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) { 209 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%s: %0.12e ", field_names[i], monitor_ctx->sum_values[i])); 210 } 211 PetscScalar sum_rates = 0; 212 for (PetscInt i = 1; i < monitor_ctx->num_comps - 1; i++) sum_rates += monitor_ctx->sum_values[i]; 213 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "SumOfRates: %0.12e ", sum_rates)); 214 PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n")); 215 } 216 PetscFunctionReturn(PETSC_SUCCESS); 217 } 218