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