xref: /libCEED/examples/fluids/problems/advection.c (revision 3e2dab0a8b85af878be609cf4667a5a492ef431d)
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 
19 PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
20   WindType             wind_type;
21   BubbleType           bubble_type;
22   BubbleContinuityType bubble_continuity_type;
23   StabilizationType    stab;
24   SetupContextAdv      setup_context;
25   User                 user = *(User *)ctx;
26   MPI_Comm             comm = user->comm;
27   Ceed                 ceed = user->ceed;
28   PetscBool            implicit;
29   AdvectionContext     advection_ctx;
30   CeedQFunctionContext advection_context;
31 
32   PetscFunctionBeginUser;
33   PetscCall(PetscCalloc1(1, &setup_context));
34   PetscCall(PetscCalloc1(1, &advection_ctx));
35 
36   // ------------------------------------------------------
37   //               SET UP ADVECTION
38   // ------------------------------------------------------
39   problem->dim                               = 3;
40   problem->q_data_size_vol                   = 10;
41   problem->q_data_size_sur                   = 10;
42   problem->setup_vol.qfunction               = Setup;
43   problem->setup_vol.qfunction_loc           = Setup_loc;
44   problem->setup_sur.qfunction               = SetupBoundary;
45   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
46   problem->ics.qfunction                     = ICsAdvection;
47   problem->ics.qfunction_loc                 = ICsAdvection_loc;
48   problem->apply_vol_rhs.qfunction           = Advection;
49   problem->apply_vol_rhs.qfunction_loc       = Advection_loc;
50   problem->apply_vol_ifunction.qfunction     = IFunction_Advection;
51   problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc;
52   problem->apply_inflow.qfunction            = Advection_InOutFlow;
53   problem->apply_inflow.qfunction_loc        = Advection_InOutFlow_loc;
54   problem->non_zero_time                     = PETSC_FALSE;
55   problem->print_info                        = PRINT_ADVECTION;
56 
57   // ------------------------------------------------------
58   //             Create the libCEED context
59   // ------------------------------------------------------
60   CeedScalar rc          = 1000.;       // m (Radius of bubble)
61   CeedScalar CtauS       = 0.;          // dimensionless
62   CeedScalar strong_form = 0.;          // [0,1]
63   CeedScalar E_wind      = 1.e6;        // J
64   PetscReal  wind[3]     = {1., 0, 0};  // m/s
65   PetscReal  domain_min[3], domain_max[3], domain_size[3];
66   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
67   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
68 
69   // ------------------------------------------------------
70   //             Create the PETSc context
71   // ------------------------------------------------------
72   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
73   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
74   PetscScalar second   = 1e-2;  // 1 second in scaled time units
75   PetscScalar Joule;
76 
77   // ------------------------------------------------------
78   //              Command line Options
79   // ------------------------------------------------------
80   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
81   // -- Physics
82   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
83   PetscBool translation;
84   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type,
85                              &translation));
86   PetscInt  n = problem->dim;
87   PetscBool user_wind;
88   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
89   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
90   PetscCall(
91       PetscOptionsScalar("-strong_form", "Strong (1) or weak/integrated by parts (0) advection residual", NULL, strong_form, &strong_form, NULL));
92   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
93   PetscCall(PetscOptionsEnum("-bubble_type", "Sphere (3D) or cylinder (2D)", NULL, BubbleTypes, (PetscEnum)(bubble_type = BUBBLE_SPHERE),
94                              (PetscEnum *)&bubble_type, NULL));
95   PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes,
96                              (PetscEnum)(bubble_continuity_type = BUBBLE_CONTINUITY_SMOOTH), (PetscEnum *)&bubble_continuity_type, NULL));
97   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
98   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
99 
100   // -- Units
101   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
102   meter = fabs(meter);
103   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
104   kilogram = fabs(kilogram);
105   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
106   second = fabs(second);
107 
108   // -- Warnings
109   if (wind_type == WIND_ROTATION && user_wind) {
110     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
111   }
112   if (wind_type == WIND_TRANSLATION && bubble_type == BUBBLE_CYLINDER && wind[2] != 0.) {
113     wind[2] = 0;
114     PetscCall(PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -bubble_type cylinder\n"));
115   }
116   if (stab == STAB_NONE && CtauS != 0) {
117     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
118   }
119   if (stab == STAB_SUPG && !implicit) {
120     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
121   }
122 
123   PetscOptionsEnd();
124 
125   // ------------------------------------------------------
126   //           Set up the PETSc context
127   // ------------------------------------------------------
128   // -- Define derived units
129   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
130 
131   user->units->meter    = meter;
132   user->units->kilogram = kilogram;
133   user->units->second   = second;
134   user->units->Joule    = Joule;
135 
136   // ------------------------------------------------------
137   //           Set up the libCEED context
138   // ------------------------------------------------------
139   // -- Scale variables to desired units
140   E_wind *= Joule;
141   rc = fabs(rc) * meter;
142   for (PetscInt i = 0; i < 3; i++) {
143     wind[i] *= (meter / second);
144     domain_size[i] *= meter;
145   }
146   problem->dm_scale = meter;
147 
148   // -- Setup Context
149   setup_context->rc                     = rc;
150   setup_context->lx                     = domain_size[0];
151   setup_context->ly                     = domain_size[1];
152   setup_context->lz                     = domain_size[2];
153   setup_context->wind[0]                = wind[0];
154   setup_context->wind[1]                = wind[1];
155   setup_context->wind[2]                = wind[2];
156   setup_context->wind_type              = wind_type;
157   setup_context->bubble_type            = bubble_type;
158   setup_context->bubble_continuity_type = bubble_continuity_type;
159   setup_context->time                   = 0;
160 
161   // -- QFunction Context
162   //  if passed correctly
163   user->phys->implicit         = implicit;
164   advection_ctx->CtauS         = CtauS;
165   advection_ctx->E_wind        = E_wind;
166   advection_ctx->implicit      = implicit;
167   advection_ctx->strong_form   = strong_form;
168   advection_ctx->stabilization = stab;
169 
170   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
171   PetscCallCeed(ceed,
172                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
173   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
174 
175   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context));
176   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
177   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc));
178   problem->apply_vol_rhs.qfunction_context = advection_context;
179   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context));
180   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context));
181   PetscFunctionReturn(PETSC_SUCCESS);
182 }
183 
184 PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx) {
185   MPI_Comm         comm = user->comm;
186   Ceed             ceed = user->ceed;
187   SetupContextAdv  setup_ctx;
188   AdvectionContext advection_ctx;
189 
190   PetscFunctionBeginUser;
191   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx));
192   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx));
193   PetscCall(PetscPrintf(comm,
194                         "  Problem:\n"
195                         "    Problem Name                       : %s\n"
196                         "    Stabilization                      : %s\n"
197                         "    Bubble Type                        : %s (%" CeedInt_FMT "D)\n"
198                         "    Bubble Continuity                  : %s\n"
199                         "    Wind Type                          : %s\n",
200                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], BubbleTypes[setup_ctx->bubble_type],
201                         setup_ctx->bubble_type == BUBBLE_SPHERE ? 3 : 2, BubbleContinuityTypes[setup_ctx->bubble_continuity_type],
202                         WindTypes[setup_ctx->wind_type]));
203 
204   if (setup_ctx->wind_type == WIND_TRANSLATION) {
205     PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2]));
206   }
207   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx));
208   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx));
209   PetscFunctionReturn(PETSC_SUCCESS);
210 }
211