xref: /libCEED/examples/fluids/problems/advection.c (revision 99e754f07c08eea4e6609a33ed68aeb9dcf08b08)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Utility functions for setting up ADVECTION
10 
11 #include "../qfunctions/advection.h"
12 
13 #include <ceed.h>
14 #include <petscdm.h>
15 
16 #include "../navierstokes.h"
17 #include "../qfunctions/setupgeo.h"
18 #include "../qfunctions/setupgeo2d.h"
19 
20 PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
21   WindType             wind_type;
22   AdvectionICType      advectionic_type;
23   BubbleContinuityType bubble_continuity_type;
24   StabilizationType    stab;
25   StabilizationTauType stab_tau;
26   SetupContextAdv      setup_context;
27   User                 user = *(User *)ctx;
28   MPI_Comm             comm = user->comm;
29   Ceed                 ceed = user->ceed;
30   PetscBool            implicit;
31   AdvectionContext     advection_ctx;
32   CeedQFunctionContext advection_context;
33   PetscInt             dim;
34 
35   PetscFunctionBeginUser;
36   PetscCall(PetscCalloc1(1, &setup_context));
37   PetscCall(PetscCalloc1(1, &advection_ctx));
38   PetscCall(DMGetDimension(dm, &dim));
39 
40   // ------------------------------------------------------
41   //               SET UP ADVECTION
42   // ------------------------------------------------------
43   switch (dim) {
44     case 2:
45       problem->dim                               = 2;
46       problem->q_data_size_vol                   = 5;
47       problem->q_data_size_sur                   = 3;
48       problem->setup_vol.qfunction               = Setup2d;
49       problem->setup_vol.qfunction_loc           = Setup2d_loc;
50       problem->setup_sur.qfunction               = SetupBoundary2d;
51       problem->setup_sur.qfunction_loc           = SetupBoundary2d_loc;
52       problem->ics.qfunction                     = ICsAdvection2d;
53       problem->ics.qfunction_loc                 = ICsAdvection2d_loc;
54       problem->apply_vol_rhs.qfunction           = RHS_Advection2d;
55       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection2d_loc;
56       problem->apply_vol_ifunction.qfunction     = IFunction_Advection2d;
57       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc;
58       problem->apply_inflow.qfunction            = Advection2d_InOutFlow;
59       problem->apply_inflow.qfunction_loc        = Advection2d_InOutFlow_loc;
60       problem->non_zero_time                     = PETSC_TRUE;
61       problem->print_info                        = PRINT_ADVECTION;
62       break;
63     case 3:
64       problem->dim                               = 3;
65       problem->q_data_size_vol                   = 10;
66       problem->q_data_size_sur                   = 10;
67       problem->setup_vol.qfunction               = Setup;
68       problem->setup_vol.qfunction_loc           = Setup_loc;
69       problem->setup_sur.qfunction               = SetupBoundary;
70       problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
71       problem->ics.qfunction                     = ICsAdvection;
72       problem->ics.qfunction_loc                 = ICsAdvection_loc;
73       problem->apply_vol_rhs.qfunction           = RHS_Advection;
74       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection_loc;
75       problem->apply_vol_ifunction.qfunction     = IFunction_Advection;
76       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc;
77       problem->apply_inflow.qfunction            = Advection_InOutFlow;
78       problem->apply_inflow.qfunction_loc        = Advection_InOutFlow_loc;
79       problem->non_zero_time                     = PETSC_FALSE;
80       problem->print_info                        = PRINT_ADVECTION;
81       break;
82   }
83 
84   // ------------------------------------------------------
85   //             Create the libCEED context
86   // ------------------------------------------------------
87   CeedScalar rc          = 1000.;  // m (Radius of bubble)
88   CeedScalar CtauS       = 0.;     // dimensionless
89   PetscBool  strong_form = PETSC_FALSE;
90   CeedScalar E_wind      = 1.e6;  // J
91   CeedScalar Ctau_a      = PetscPowScalarInt(user->app_ctx->degree, 2);
92   CeedScalar Ctau_t      = 0.;
93   PetscReal  wind[3]     = {1., 0, 0};  // m/s
94   PetscReal  domain_min[3], domain_max[3], domain_size[3];
95   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
96   for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
97 
98   // ------------------------------------------------------
99   //             Create the PETSc context
100   // ------------------------------------------------------
101   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
102   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
103   PetscScalar second   = 1e-2;  // 1 second in scaled time units
104   PetscScalar Joule;
105 
106   // ------------------------------------------------------
107   //              Command line Options
108   // ------------------------------------------------------
109   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
110   // -- Physics
111   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
112   PetscBool translation;
113   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type,
114                              &translation));
115   PetscInt  n = problem->dim;
116   PetscBool user_wind;
117   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
118   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
119   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
120                              NULL));
121   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
122   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes,
123                              (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
124   bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE;
125   PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type,
126                              (PetscEnum *)&bubble_continuity_type, NULL));
127   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
128   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
129                              (PetscEnum *)&stab_tau, NULL));
130   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
131   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL));
132   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
133 
134   // -- Units
135   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
136   meter = fabs(meter);
137   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
138   kilogram = fabs(kilogram);
139   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
140   second = fabs(second);
141 
142   // -- Warnings
143   if (wind_type == WIND_ROTATION && user_wind) {
144     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
145   }
146   if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) {
147     wind[2] = 0;
148     PetscCall(
149         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
150   }
151   if (stab == STAB_NONE && CtauS != 0) {
152     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
153   }
154   if (stab == STAB_SUPG && !implicit) {
155     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
156   }
157 
158   PetscOptionsEnd();
159 
160   // ------------------------------------------------------
161   //           Set up the PETSc context
162   // ------------------------------------------------------
163   // -- Define derived units
164   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
165 
166   user->units->meter    = meter;
167   user->units->kilogram = kilogram;
168   user->units->second   = second;
169   user->units->Joule    = Joule;
170 
171   // ------------------------------------------------------
172   //           Set up the libCEED context
173   // ------------------------------------------------------
174   // -- Scale variables to desired units
175   E_wind *= Joule;
176   rc = fabs(rc) * meter;
177   for (PetscInt i = 0; i < problem->dim; i++) {
178     wind[i] *= (meter / second);
179     domain_size[i] *= meter;
180   }
181   problem->dm_scale = meter;
182 
183   // -- Setup Context
184   setup_context->rc                     = rc;
185   setup_context->lx                     = domain_size[0];
186   setup_context->ly                     = domain_size[1];
187   setup_context->lz                     = problem->dim == 3 ? domain_size[2] : 0.;
188   setup_context->wind[0]                = wind[0];
189   setup_context->wind[1]                = wind[1];
190   setup_context->wind[2]                = problem->dim == 3 ? wind[2] : 0.;
191   setup_context->wind_type              = wind_type;
192   setup_context->initial_condition_type = advectionic_type;
193   setup_context->bubble_continuity_type = bubble_continuity_type;
194   setup_context->time                   = 0;
195 
196   // -- QFunction Context
197   user->phys->implicit             = implicit;
198   advection_ctx->CtauS             = CtauS;
199   advection_ctx->E_wind            = E_wind;
200   advection_ctx->implicit          = implicit;
201   advection_ctx->strong_form       = strong_form;
202   advection_ctx->stabilization     = stab;
203   advection_ctx->stabilization_tau = stab_tau;
204   advection_ctx->Ctau_a            = Ctau_a;
205   advection_ctx->Ctau_t            = Ctau_t;
206 
207   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
208   PetscCallCeed(ceed,
209                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
210   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
211 
212   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context));
213   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
214   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc));
215   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
216                                                          "Size of timestep, delta t"));
217   problem->apply_vol_rhs.qfunction_context = advection_context;
218   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context));
219   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context));
220   PetscFunctionReturn(PETSC_SUCCESS);
221 }
222 
223 PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx) {
224   MPI_Comm         comm = user->comm;
225   Ceed             ceed = user->ceed;
226   SetupContextAdv  setup_ctx;
227   AdvectionContext advection_ctx;
228 
229   PetscFunctionBeginUser;
230   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx));
231   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx));
232   PetscCall(PetscPrintf(comm,
233                         "  Problem:\n"
234                         "    Problem Name                       : %s\n"
235                         "    Stabilization                      : %s\n"
236                         "    Initial Condition Type             : %s\n"
237                         "    Bubble Continuity                  : %s\n"
238                         "    Wind Type                          : %s\n",
239                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type],
240                         BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type]));
241 
242   if (setup_ctx->wind_type == WIND_TRANSLATION) {
243     switch (problem->dim) {
244       case 2:
245         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1]));
246         break;
247       case 3:
248         PetscCall(
249             PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2]));
250         break;
251     }
252   }
253   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx));
254   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx));
255   PetscFunctionReturn(PETSC_SUCCESS);
256 }
257