125125139SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
225125139SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
325125139SJames Wright #include "../qfunctions/monitor_totalkineticenergy.h"
4a5677b81SMohi #include <honee.h>
525125139SJames Wright #include <navierstokes.h>
625125139SJames Wright #include <time.h>
725125139SJames Wright
825125139SJames Wright typedef struct {
925125139SJames Wright OperatorApplyContext op_monitor_ctx;
1025125139SJames Wright Vec values;
1125125139SJames Wright PetscScalar *sum_values;
1225125139SJames Wright PetscScalar volume;
1325125139SJames Wright PetscInt num_comps;
142a51b432SJames Wright PetscInt tab_level;
1525125139SJames Wright PetscBool is_header_written;
16600dae7dSJames Wright } *MonitorTotalKE;
1725125139SJames Wright
MonitorTotalKEDestroy(MonitorTotalKE * monitor_ctx)185206a5a0SJames Wright static PetscErrorCode MonitorTotalKEDestroy(MonitorTotalKE *monitor_ctx) {
195206a5a0SJames Wright MonitorTotalKE monitor_ctx_ = *monitor_ctx;
2025125139SJames Wright
2125125139SJames Wright PetscFunctionBeginUser;
225206a5a0SJames Wright PetscCall(OperatorApplyContextDestroy(monitor_ctx_->op_monitor_ctx));
235206a5a0SJames Wright PetscCall(VecDestroy(&monitor_ctx_->values));
245206a5a0SJames Wright PetscCall(PetscFree(monitor_ctx_->sum_values));
255206a5a0SJames Wright PetscCall(PetscFree(monitor_ctx_));
2625125139SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
2725125139SJames Wright }
2825125139SJames Wright
2925125139SJames Wright /** @brief Create CeedElemRestriction for collocated data in component-major order.
3025125139SJames Wright a. Sets the strides of the restriction to component-major order
3125125139SJames Wright Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
3225125139SJames Wright */
CreateElemRestrColloc_CompMajor(Ceed ceed,CeedInt num_comp,CeedBasis basis,CeedElemRestriction elem_restr_base,CeedElemRestriction * elem_restr_collocated)3325125139SJames Wright static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
3425125139SJames Wright CeedElemRestriction *elem_restr_collocated) {
3525125139SJames Wright CeedInt num_elem_qpts, loc_num_elem;
3625125139SJames Wright
3725125139SJames Wright PetscFunctionBeginUser;
3825125139SJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
3925125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
4025125139SJames Wright
4125125139SJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
4225125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
4325125139SJames Wright elem_restr_collocated));
4425125139SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
4525125139SJames Wright }
4625125139SJames Wright
SetupMontiorTotalKineticEnergy(TS ts,PetscViewerAndFormat * ctx)4725125139SJames Wright PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx) {
4825125139SJames Wright Honee honee;
4925125139SJames Wright Ceed ceed;
5025125139SJames Wright CeedQFunction qf_monitor;
5125125139SJames Wright CeedOperator op_monitor;
5225125139SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_totalke, elem_restr_q;
5325125139SJames Wright CeedBasis basis_q;
5425125139SJames Wright CeedVector q_data;
5525125139SJames Wright CeedInt num_comp_q, q_data_size;
56e3db12f8SJames Wright PetscInt num_comp_totalke = 5, dim, tab_level;
57600dae7dSJames Wright MonitorTotalKE monitor_ctx;
5825125139SJames Wright CeedQFunctionContext newt_qfctx;
5925125139SJames Wright PetscBool is_ascii;
6025125139SJames Wright
6125125139SJames Wright PetscFunctionBeginUser;
6225125139SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii));
6325125139SJames Wright PetscCheck(is_ascii, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only supports ASCII viewers");
6425125139SJames Wright PetscCall(TSGetApplicationContext(ts, &honee));
6525125139SJames Wright PetscCall(DMGetDimension(honee->dm, &dim));
6625125139SJames Wright ceed = honee->ceed;
6725125139SJames Wright
6825125139SJames Wright PetscCall(PetscNew(&monitor_ctx));
6925125139SJames Wright PetscCall(HoneeCalculateDomainSize(honee, &monitor_ctx->volume));
7025125139SJames Wright { // Get restriction and basis from the RHS function
716f5e6828SJames Wright CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
7225125139SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op
7325125139SJames Wright
74da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(main_op, &sub_ops));
7525125139SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newt_qfctx));
7625125139SJames Wright }
77*9018c49aSJames Wright PetscCall(QDataGet(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
78cf8f12c8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
79cf8f12c8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
8025125139SJames Wright PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_totalke, basis_q, elem_restr_q, &elem_restr_totalke));
81cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
8225125139SJames Wright
8325125139SJames Wright switch (honee->phys->state_var) {
8425125139SJames Wright case STATEVAR_PRIMITIVE:
8525125139SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Prim, MonitorTotalKineticEnergy_Prim_loc, &qf_monitor));
8625125139SJames Wright break;
8725125139SJames Wright case STATEVAR_CONSERVATIVE:
8825125139SJames Wright PetscCallCeed(ceed,
8925125139SJames Wright CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Conserv, MonitorTotalKineticEnergy_Conserv_loc, &qf_monitor));
9025125139SJames Wright break;
9125125139SJames Wright case STATEVAR_ENTROPY:
9225125139SJames Wright PetscCallCeed(ceed,
9325125139SJames Wright CeedQFunctionCreateInterior(ceed, 1, MonitorTotalKineticEnergy_Entropy, MonitorTotalKineticEnergy_Entropy_loc, &qf_monitor));
9425125139SJames Wright break;
9525125139SJames Wright }
9625125139SJames Wright
9725125139SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_monitor, newt_qfctx));
9825125139SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q", num_comp_q, CEED_EVAL_INTERP));
9925125139SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
10025125139SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q_data", q_data_size, CEED_EVAL_NONE));
10125125139SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_monitor, "v", num_comp_totalke, CEED_EVAL_NONE));
10225125139SJames Wright
10325125139SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_monitor, NULL, NULL, &op_monitor));
10425125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
10525125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
10625125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data));
10725125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "v", elem_restr_totalke, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
10825125139SJames Wright
10925125139SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, NULL, ceed, op_monitor, NULL, NULL, NULL, NULL, &monitor_ctx->op_monitor_ctx));
11025125139SJames Wright
11125125139SJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_monitor, DMReturnVecType(honee->dm), PETSC_COMM_SELF, NULL, &monitor_ctx->values));
11225125139SJames Wright PetscCall(VecSetBlockSize(monitor_ctx->values, num_comp_totalke));
11325125139SJames Wright monitor_ctx->num_comps = num_comp_totalke;
11425125139SJames Wright PetscCall(PetscMalloc1(num_comp_totalke, &monitor_ctx->sum_values));
1152a51b432SJames Wright PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tab_level));
1162a51b432SJames Wright monitor_ctx->tab_level = tab_level + 1;
11725125139SJames Wright
11825125139SJames Wright ctx->data = monitor_ctx;
1195206a5a0SJames Wright ctx->data_destroy = (PetscCtxDestroyFn *)MonitorTotalKEDestroy;
12025125139SJames Wright
12125125139SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newt_qfctx));
12225125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
12325125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
12425125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_totalke));
12525125139SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
12625125139SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
12725125139SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_monitor));
12825125139SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_monitor));
12925125139SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
13025125139SJames Wright }
13125125139SJames Wright
TSMonitor_TotalKineticEnergy(TS ts,PetscInt step,PetscReal solution_time,Vec Q,PetscViewerAndFormat * ctx)13225125139SJames Wright PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
133600dae7dSJames Wright MonitorTotalKE monitor_ctx = (MonitorTotalKE)ctx->data;
13425125139SJames Wright Honee honee;
13525125139SJames Wright MPI_Comm comm;
13625125139SJames Wright PetscMPIInt rank;
13725125139SJames Wright TSConvergedReason reason;
13825125139SJames Wright static const char *field_names[] = {"TotalKineticEnergy", "StrainDissipation", "DivergenceDissipation", "VolumeExpansion", "MuVorticity2"};
13925125139SJames Wright
14025125139SJames Wright PetscFunctionBeginUser;
14125125139SJames Wright PetscCall(TSGetConvergedReason(ts, &reason));
14225125139SJames Wright if (!(step % ctx->view_interval == 0 || reason != TS_CONVERGED_ITERATING)) PetscFunctionReturn(PETSC_SUCCESS);
14325125139SJames Wright PetscCall(TSGetApplicationContext(ts, &honee));
14425125139SJames Wright PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
14525125139SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank));
14625125139SJames Wright
14725125139SJames Wright PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time));
14825125139SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc));
14925125139SJames Wright
15025125139SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(honee->Q_loc, monitor_ctx->values, monitor_ctx->op_monitor_ctx));
15125125139SJames Wright
15225125139SJames Wright for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) {
15325125139SJames Wright PetscCall(VecStrideSum(monitor_ctx->values, i, &monitor_ctx->sum_values[i]));
15425125139SJames Wright }
15525125139SJames Wright
15625125139SJames Wright if (rank == 0) PetscCallMPI(MPI_Reduce(MPI_IN_PLACE, monitor_ctx->sum_values, monitor_ctx->num_comps, MPIU_SCALAR, MPI_SUM, 0, comm));
15725125139SJames Wright else PetscCallMPI(MPI_Reduce(monitor_ctx->sum_values, monitor_ctx->sum_values, monitor_ctx->num_comps, MPIU_SCALAR, MPI_SUM, 0, comm));
15825125139SJames Wright
15925125139SJames Wright if (rank == 0)
16025125139SJames Wright for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) monitor_ctx->sum_values[i] /= monitor_ctx->volume;
16125125139SJames Wright
16225125139SJames Wright if (ctx->format == PETSC_VIEWER_ASCII_CSV) {
16325125139SJames Wright if (!monitor_ctx->is_header_written) {
16425125139SJames Wright char buf[PETSC_MAX_PATH_LEN];
16525125139SJames Wright const char *buf_const;
166a5677b81SMohi int major, minor, patch;
16725125139SJames Wright time_t t = time(NULL);
16825125139SJames Wright
169a5677b81SMohi PetscCall(HoneeGetVersion(&major, &minor, &patch, NULL));
170a5677b81SMohi PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# honee_version: %d.%d.%d\n", major, minor, patch));
171a5677b81SMohi PetscCall(HoneeGetGitVersion(&buf_const));
172a5677b81SMohi PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# honee_git_commit: %s\n", buf_const));
17325125139SJames Wright PetscCall(PetscGetVersion(buf, sizeof(buf)));
17425125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_version: %s\n", buf));
17525125139SJames Wright PetscCall(PetscGetPetscDir(&buf_const));
17625125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_directory: %s\n", buf_const));
17725125139SJames Wright PetscCall(PetscGetArchType(buf, sizeof(buf)));
17825125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_arch: %s\n", buf));
17914bd2a07SJames Wright PetscCheck(strftime(buf, sizeof(buf), "%FT%T%z", localtime(&t)), comm, PETSC_ERR_SYS, "strftime() call failed");
18025125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date: %s\n", buf));
18114bd2a07SJames Wright PetscCheck(strftime(buf, sizeof(buf), "%Z", localtime(&t)), comm, PETSC_ERR_SYS, "strftime() call failed");
18225125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date_timezone: %s\n", buf));
18325125139SJames Wright PetscCall(PetscGetUserName(buf, sizeof(buf)));
18425125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# username: %s\n", buf));
18525125139SJames Wright PetscCall(PetscGetHostName(buf, sizeof(buf)));
18625125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# hostname: %s\n", buf));
18725125139SJames Wright PetscCall(PetscGetWorkingDirectory(buf, sizeof(buf)));
18825125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# working_directory: %s\n", buf));
18925125139SJames Wright PetscCall(PetscViewerFileGetName(ctx->viewer, &buf_const));
19025125139SJames Wright PetscCall(PetscGetFullPath(buf_const, buf, sizeof(buf)));
19125125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# original_file_path: %s\n", buf));
19225125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "#\n"));
19325125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "step,time,"));
19425125139SJames Wright for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) {
19525125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%s", field_names[i]));
19625125139SJames Wright if (i < monitor_ctx->num_comps - 1) PetscCall(PetscViewerASCIIPrintf(ctx->viewer, ","));
19725125139SJames Wright }
19825125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n"));
19925125139SJames Wright monitor_ctx->is_header_written = PETSC_TRUE;
20025125139SJames Wright }
20125125139SJames Wright
20225125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%" PetscInt_FMT ",%0.17e,", step, solution_time));
20325125139SJames Wright for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) {
20425125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%0.17e", monitor_ctx->sum_values[i]));
20525125139SJames Wright if (i < monitor_ctx->num_comps - 1) PetscCall(PetscViewerASCIIPrintf(ctx->viewer, ","));
20625125139SJames Wright }
20725125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n"));
20825125139SJames Wright } else {
2092a51b432SJames Wright PetscCall(PetscViewerASCIIAddTab(ctx->viewer, monitor_ctx->tab_level));
21025125139SJames Wright for (PetscInt i = 0; i < monitor_ctx->num_comps; i++) {
21125125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%s: %0.12e ", field_names[i], monitor_ctx->sum_values[i]));
2122a51b432SJames Wright if (i == 0) PetscCall(PetscViewerASCIIUseTabs(ctx->viewer, PETSC_FALSE));
21325125139SJames Wright }
21425125139SJames Wright PetscScalar sum_rates = 0;
21525125139SJames Wright for (PetscInt i = 1; i < monitor_ctx->num_comps - 1; i++) sum_rates += monitor_ctx->sum_values[i];
21625125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "SumOfRates: %0.12e ", sum_rates));
21725125139SJames Wright PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n"));
2182a51b432SJames Wright PetscCall(PetscViewerASCIIUseTabs(ctx->viewer, PETSC_TRUE));
2192a51b432SJames Wright PetscCall(PetscViewerASCIISubtractTab(ctx->viewer, monitor_ctx->tab_level));
22025125139SJames Wright }
22125125139SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
22225125139SJames Wright }
223