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