xref: /honee/src/monitor_cfl.c (revision eea329a2b230f0fd8b62e58c0441b5eed4fc5ca5)
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_cfl.h"
5 #include <navierstokes.h>
6 
7 typedef struct {
8   OperatorApplyContext op_monitor_ctx;
9   Vec                  values;
10   PetscInt             num_comps;
11   PetscInt             tab_level;
12   PetscBool            is_header_written;
13 } *MonitorCflCtx;
14 
15 PetscErrorCode MonitorCflCtxDestroy(void **ctx) {
16   MonitorCflCtx monitor_ctx = *(MonitorCflCtx *)ctx;
17 
18   PetscFunctionBeginUser;
19   PetscCall(OperatorApplyContextDestroy(monitor_ctx->op_monitor_ctx));
20   PetscCall(VecDestroy(&monitor_ctx->values));
21   PetscCall(PetscFree(monitor_ctx));
22   PetscFunctionReturn(PETSC_SUCCESS);
23 }
24 
25 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx) {
26   Honee                honee;
27   Ceed                 ceed;
28   CeedQFunction        qf_monitor = NULL;
29   CeedOperator         op_monitor;
30   CeedElemRestriction  elem_restr_qd, elem_restr_cfl, elem_restr_q;
31   CeedBasis            basis_q;
32   CeedVector           q_data;
33   CeedInt              num_comp_q, q_data_size;
34   DMLabel              domain_label = NULL;
35   PetscInt             label_value = 0, num_comp_cfl = 1, dim, tab_level;
36   MonitorCflCtx        monitor_ctx;
37   CeedQFunctionContext newt_qfctx;
38   PetscBool            is_ascii;
39 
40   PetscFunctionBeginUser;
41   PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii));
42   PetscCheck(is_ascii, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Only supports ASCII viewers");
43   PetscCall(TSGetApplicationContext(ts, &honee));
44   PetscCall(DMGetDimension(honee->dm, &dim));
45   ceed = honee->ceed;
46 
47   PetscCall(PetscNew(&monitor_ctx));
48 
49   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
50   PetscCall(QDataGet(ceed, honee->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
51                      &q_data_size));
52   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, honee->dm, domain_label, label_value, 0, num_comp_cfl, &elem_restr_cfl));
53 
54   {  // Get restriction and basis from the RHS function
55     CeedOperator     *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
56     CeedOperatorField op_field;
57     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
58 
59     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(main_op, &sub_ops));
60     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
61 
62     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
63     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newt_qfctx));
64   }
65 
66   switch (dim) {
67     case 2:
68       switch (honee->phys->state_var) {
69         case STATEVAR_PRIMITIVE:
70           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_2D_Prim, MonitorCFL_2D_Prim_loc, &qf_monitor));
71           break;
72         case STATEVAR_CONSERVATIVE:
73           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_2D_Conserv, MonitorCFL_2D_Conserv_loc, &qf_monitor));
74           break;
75         case STATEVAR_ENTROPY:
76           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_2D_Entropy, MonitorCFL_2D_Entropy_loc, &qf_monitor));
77           break;
78       }
79       break;
80     case 3:
81       switch (honee->phys->state_var) {
82         case STATEVAR_PRIMITIVE:
83           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_3D_Prim, MonitorCFL_3D_Prim_loc, &qf_monitor));
84           break;
85         case STATEVAR_CONSERVATIVE:
86           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_3D_Conserv, MonitorCFL_3D_Conserv_loc, &qf_monitor));
87           break;
88         case STATEVAR_ENTROPY:
89           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MonitorCFL_3D_Entropy, MonitorCFL_3D_Entropy_loc, &qf_monitor));
90           break;
91       }
92       break;
93   }
94   PetscCheck(qf_monitor, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP,
95              "Could not create CFL monitor QFunction for dim %" PetscInt_FMT " and state variable %s", dim, StateVariables[honee->phys->state_var]);
96 
97   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_monitor, newt_qfctx));
98   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q", num_comp_q, CEED_EVAL_INTERP));
99   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_monitor, "q_data", q_data_size, CEED_EVAL_NONE));
100   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_monitor, "v", num_comp_cfl, CEED_EVAL_NONE));
101 
102   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_monitor, NULL, NULL, &op_monitor));
103   PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
104   PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data));
105   PetscCallCeed(ceed, CeedOperatorSetField(op_monitor, "v", elem_restr_cfl, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
106 
107   PetscCall(OperatorApplyContextCreate(honee->dm, NULL, ceed, op_monitor, NULL, NULL, NULL, NULL, &monitor_ctx->op_monitor_ctx));
108 
109   PetscCall(CeedOperatorCreateLocalVecs(op_monitor, DMReturnVecType(honee->dm), PETSC_COMM_SELF, NULL, &monitor_ctx->values));
110   PetscCall(VecSetBlockSize(monitor_ctx->values, num_comp_cfl));
111   monitor_ctx->num_comps = num_comp_cfl;
112   PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tab_level));
113   monitor_ctx->tab_level = tab_level + 1;
114 
115   ctx->data         = monitor_ctx;
116   ctx->data_destroy = MonitorCflCtxDestroy;
117 
118   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newt_qfctx));
119   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
120   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
121   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_cfl));
122   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
123   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
124   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_monitor));
125   PetscCallCeed(ceed, CeedOperatorDestroy(&op_monitor));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
130   MonitorCflCtx     monitor_ctx = (MonitorCflCtx)ctx->data;
131   Honee             honee;
132   MPI_Comm          comm;
133   TSConvergedReason reason;
134   PetscReal         part_minmax[2], global_minmax[2], dt;
135 
136   PetscFunctionBeginUser;
137   PetscCall(TSGetConvergedReason(ts, &reason));
138   if (!(step % ctx->view_interval == 0 || reason != TS_CONVERGED_ITERATING)) PetscFunctionReturn(PETSC_SUCCESS);
139   PetscCall(TSGetApplicationContext(ts, &honee));
140   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
141 
142   PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time));
143   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc));
144 
145   PetscCall(ApplyCeedOperatorLocalToLocal(honee->Q_loc, monitor_ctx->values, monitor_ctx->op_monitor_ctx));
146 
147   PetscCall(VecMin(monitor_ctx->values, NULL, &part_minmax[0]));
148   PetscCall(VecMax(monitor_ctx->values, NULL, &part_minmax[1]));
149   // Need to get global min/max since `values` is only the local partition
150   PetscCall(PetscGlobalMinMaxReal(comm, part_minmax, global_minmax));
151   PetscCall(TSGetTimeStep(ts, &dt));
152   global_minmax[0] *= dt;
153   global_minmax[1] *= dt;
154 
155   if (ctx->format == PETSC_VIEWER_ASCII_CSV) {
156     if (!monitor_ctx->is_header_written) {
157       char        buf[PETSC_MAX_PATH_LEN];
158       const char *buf_const;
159       time_t      t = time(NULL);
160 
161       PetscCall(PetscGetVersion(buf, sizeof(buf)));
162       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_version: %s\n", buf));
163       PetscCall(PetscGetPetscDir(&buf_const));
164       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_directory: %s\n", buf_const));
165       PetscCall(PetscGetArchType(buf, sizeof(buf)));
166       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# petsc_arch: %s\n", buf));
167       if (strftime(buf, sizeof(buf), "%FT%T%z", localtime(&t)) == 0) SETERRQ(comm, PETSC_ERR_SYS, "strftime() call failed");
168       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date: %s\n", buf));
169       if (strftime(buf, sizeof(buf), "%Z", localtime(&t)) == 0) SETERRQ(comm, PETSC_ERR_SYS, "strftime() call failed");
170       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# date_timezone: %s\n", buf));
171       PetscCall(PetscGetUserName(buf, sizeof(buf)));
172       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# username: %s\n", buf));
173       PetscCall(PetscGetHostName(buf, sizeof(buf)));
174       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# hostname: %s\n", buf));
175       PetscCall(PetscGetWorkingDirectory(buf, sizeof(buf)));
176       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# working_directory: %s\n", buf));
177 
178       PetscCall(PetscViewerFileGetName(ctx->viewer, &buf_const));
179       PetscCall(PetscGetFullPath(buf_const, buf, sizeof(buf)));
180       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "# original_file_path: %s\n", buf));
181       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "#\n"));
182       PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "step,time,mincfl,maxcfl\n"));
183       monitor_ctx->is_header_written = PETSC_TRUE;
184     }
185 
186     PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%" PetscInt_FMT ",%0.17e,", step, solution_time));
187     PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "%0.17e,%0.17e", global_minmax[0], global_minmax[1]));
188     PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "\n"));
189   } else {
190     PetscCall(PetscViewerASCIIAddTab(ctx->viewer, monitor_ctx->tab_level));
191     PetscCall(PetscViewerASCIIPrintf(ctx->viewer, "Min, Max CFL: %0.5g, %0.5g\n", global_minmax[0], global_minmax[1]));
192     PetscCall(PetscViewerASCIISubtractTab(ctx->viewer, monitor_ctx->tab_level));
193   }
194   PetscFunctionReturn(PETSC_SUCCESS);
195 }
196