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