1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details. 4ccaff030SJeremy L Thompson // 5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software 6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral 7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and 8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed. 9ccaff030SJeremy L Thompson // 10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office 12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for 13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including 14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early 15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative. 16ccaff030SJeremy L Thompson 17ccaff030SJeremy L Thompson // libCEED + PETSc Example: Navier-Stokes 18ccaff030SJeremy L Thompson // 19ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve a 20ccaff030SJeremy L Thompson // Navier-Stokes problem. 21ccaff030SJeremy L Thompson // 22ccaff030SJeremy L Thompson // The code is intentionally "raw", using only low-level communication 23ccaff030SJeremy L Thompson // primitives. 24ccaff030SJeremy L Thompson // 25ccaff030SJeremy L Thompson // Build with: 26ccaff030SJeremy L Thompson // 27ccaff030SJeremy L Thompson // make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes 28ccaff030SJeremy L Thompson // 29ccaff030SJeremy L Thompson // Sample runs: 30ccaff030SJeremy L Thompson // 31ff6701fcSJed Brown // ./navierstokes -ceed /cpu/self -problem density_current -degree 1 32ff6701fcSJed Brown // ./navierstokes -ceed /gpu/occa -problem advection -degree 1 33ccaff030SJeremy L Thompson // 348ca749c7Svaleriabarra //TESTARGS -ceed {ceed_resource} -test explicit -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 358ca749c7Svaleriabarra //TESTARGS -ceed {ceed_resource} -test implicit_stab_none -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha 368ca749c7Svaleriabarra //TESTARGS -ceed {ceed_resource} -test implicit_stab_supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -stab supg 37ccaff030SJeremy L Thompson 38ccaff030SJeremy L Thompson /// @file 39ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc 40ccaff030SJeremy L Thompson 41ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 42ccaff030SJeremy L Thompson 43ccaff030SJeremy L Thompson #include <petscts.h> 44ccaff030SJeremy L Thompson #include <petscdmplex.h> 45ccaff030SJeremy L Thompson #include <ceed.h> 46ccaff030SJeremy L Thompson #include <stdbool.h> 47ccaff030SJeremy L Thompson #include <petscsys.h> 48ccaff030SJeremy L Thompson #include "common.h" 49ccaff030SJeremy L Thompson #include "advection.h" 50ccaff030SJeremy L Thompson #include "advection2d.h" 51ccaff030SJeremy L Thompson #include "densitycurrent.h" 52ccaff030SJeremy L Thompson 530261bf83Svaleriabarra // MemType Options 540261bf83Svaleriabarra static const char *const memTypes[] = { 550261bf83Svaleriabarra "host", 560261bf83Svaleriabarra "device", 570261bf83Svaleriabarra "memType", "CEED_MEM_", NULL 580261bf83Svaleriabarra }; 590261bf83Svaleriabarra 60ccaff030SJeremy L Thompson // Problem Options 61ccaff030SJeremy L Thompson typedef enum { 62ccaff030SJeremy L Thompson NS_DENSITY_CURRENT = 0, 63ccaff030SJeremy L Thompson NS_ADVECTION = 1, 64ccaff030SJeremy L Thompson NS_ADVECTION2D = 2, 65ccaff030SJeremy L Thompson } problemType; 66ccaff030SJeremy L Thompson static const char *const problemTypes[] = { 67ccaff030SJeremy L Thompson "density_current", 68ccaff030SJeremy L Thompson "advection", 69ccaff030SJeremy L Thompson "advection2d", 702f9599b8Svaleriabarra "problemType", "NS_", NULL 71ccaff030SJeremy L Thompson }; 72ccaff030SJeremy L Thompson 73ccaff030SJeremy L Thompson typedef enum { 74ccaff030SJeremy L Thompson STAB_NONE = 0, 75ccaff030SJeremy L Thompson STAB_SU = 1, // Streamline Upwind 76ccaff030SJeremy L Thompson STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 77ccaff030SJeremy L Thompson } StabilizationType; 78ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = { 79*6f55dfd5Svaleriabarra "none", 80ccaff030SJeremy L Thompson "SU", 81ccaff030SJeremy L Thompson "SUPG", 82ccaff030SJeremy L Thompson "StabilizationType", "STAB_", NULL 83ccaff030SJeremy L Thompson }; 84ccaff030SJeremy L Thompson 858ca749c7Svaleriabarra // Test Options 868ca749c7Svaleriabarra typedef enum { 878ca749c7Svaleriabarra TEST_NONE = 0, // Non test mode 888ca749c7Svaleriabarra TEST_EXPLICIT = 1, // Explicit test 898ca749c7Svaleriabarra TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab 908ca749c7Svaleriabarra TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab 918ca749c7Svaleriabarra } testType; 928ca749c7Svaleriabarra static const char *const testTypes[] = { 938ca749c7Svaleriabarra "none", 948ca749c7Svaleriabarra "explicit", 958ca749c7Svaleriabarra "implicit_stab_none", 968ca749c7Svaleriabarra "implicit_stab_supg", 978ca749c7Svaleriabarra "testType", "TEST_", NULL 988ca749c7Svaleriabarra }; 998ca749c7Svaleriabarra 1008ca749c7Svaleriabarra // Tests specific data 1018ca749c7Svaleriabarra typedef struct { 1028ca749c7Svaleriabarra PetscScalar testtol; 1038ca749c7Svaleriabarra const char *filepath; 1048ca749c7Svaleriabarra } testData; 1058ca749c7Svaleriabarra 1068ca749c7Svaleriabarra testData testOptions[] = { 1078ca749c7Svaleriabarra [TEST_NONE] = { 1088ca749c7Svaleriabarra .testtol = 0., 1098ca749c7Svaleriabarra .filepath = NULL 1108ca749c7Svaleriabarra }, 1118ca749c7Svaleriabarra [TEST_EXPLICIT] = { 1128ca749c7Svaleriabarra .testtol = 1E-5, 1138ca749c7Svaleriabarra .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin" 1148ca749c7Svaleriabarra }, 1158ca749c7Svaleriabarra [TEST_IMPLICIT_STAB_NONE] = { 1168ca749c7Svaleriabarra .testtol = 5E-4, 1178ca749c7Svaleriabarra .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin" 1188ca749c7Svaleriabarra }, 1198ca749c7Svaleriabarra [TEST_IMPLICIT_STAB_SUPG] = { 1208ca749c7Svaleriabarra .testtol = 5E-4, 1218ca749c7Svaleriabarra .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin" 1228ca749c7Svaleriabarra } 1238ca749c7Svaleriabarra }; 1248ca749c7Svaleriabarra 125ccaff030SJeremy L Thompson // Problem specific data 126ccaff030SJeremy L Thompson typedef struct { 127ccaff030SJeremy L Thompson CeedInt dim, qdatasize; 128ccaff030SJeremy L Thompson CeedQFunctionUser setup, ics, apply_rhs, apply_ifunction; 129ccaff030SJeremy L Thompson PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 130ccaff030SJeremy L Thompson PetscScalar[], void *); 131ccaff030SJeremy L Thompson const char *setup_loc, *ics_loc, *apply_rhs_loc, *apply_ifunction_loc; 132ccaff030SJeremy L Thompson const bool non_zero_time; 133ccaff030SJeremy L Thompson } problemData; 134ccaff030SJeremy L Thompson 135ccaff030SJeremy L Thompson problemData problemOptions[] = { 136ccaff030SJeremy L Thompson [NS_DENSITY_CURRENT] = { 137ccaff030SJeremy L Thompson .dim = 3, 138ccaff030SJeremy L Thompson .qdatasize = 10, 139ccaff030SJeremy L Thompson .setup = Setup, 140ccaff030SJeremy L Thompson .setup_loc = Setup_loc, 141ccaff030SJeremy L Thompson .ics = ICsDC, 142ccaff030SJeremy L Thompson .ics_loc = ICsDC_loc, 143ccaff030SJeremy L Thompson .apply_rhs = DC, 144ccaff030SJeremy L Thompson .apply_rhs_loc = DC_loc, 145ccaff030SJeremy L Thompson .apply_ifunction = IFunction_DC, 146ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_DC_loc, 147ccaff030SJeremy L Thompson .bc = Exact_DC, 1483e4faa5aSvaleriabarra .non_zero_time = PETSC_FALSE, 149ccaff030SJeremy L Thompson }, 150ccaff030SJeremy L Thompson [NS_ADVECTION] = { 151ccaff030SJeremy L Thompson .dim = 3, 152ccaff030SJeremy L Thompson .qdatasize = 10, 153ccaff030SJeremy L Thompson .setup = Setup, 154ccaff030SJeremy L Thompson .setup_loc = Setup_loc, 155ccaff030SJeremy L Thompson .ics = ICsAdvection, 156ccaff030SJeremy L Thompson .ics_loc = ICsAdvection_loc, 157ccaff030SJeremy L Thompson .apply_rhs = Advection, 158ccaff030SJeremy L Thompson .apply_rhs_loc = Advection_loc, 159ccaff030SJeremy L Thompson .apply_ifunction = IFunction_Advection, 160ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_Advection_loc, 161ccaff030SJeremy L Thompson .bc = Exact_Advection, 1623e4faa5aSvaleriabarra .non_zero_time = PETSC_FALSE, 163ccaff030SJeremy L Thompson }, 164ccaff030SJeremy L Thompson [NS_ADVECTION2D] = { 165ccaff030SJeremy L Thompson .dim = 2, 166ccaff030SJeremy L Thompson .qdatasize = 5, 167ccaff030SJeremy L Thompson .setup = Setup2d, 168ccaff030SJeremy L Thompson .setup_loc = Setup2d_loc, 169ccaff030SJeremy L Thompson .ics = ICsAdvection2d, 170ccaff030SJeremy L Thompson .ics_loc = ICsAdvection2d_loc, 171ccaff030SJeremy L Thompson .apply_rhs = Advection2d, 172ccaff030SJeremy L Thompson .apply_rhs_loc = Advection2d_loc, 173ccaff030SJeremy L Thompson .apply_ifunction = IFunction_Advection2d, 174ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_Advection2d_loc, 175ccaff030SJeremy L Thompson .bc = Exact_Advection2d, 1763e4faa5aSvaleriabarra .non_zero_time = PETSC_TRUE, 177ccaff030SJeremy L Thompson }, 178ccaff030SJeremy L Thompson }; 179ccaff030SJeremy L Thompson 180ccaff030SJeremy L Thompson // PETSc user data 181ccaff030SJeremy L Thompson typedef struct User_ *User; 182ccaff030SJeremy L Thompson typedef struct Units_ *Units; 183ccaff030SJeremy L Thompson 184ccaff030SJeremy L Thompson struct User_ { 185ccaff030SJeremy L Thompson MPI_Comm comm; 186ccaff030SJeremy L Thompson PetscInt outputfreq; 187ccaff030SJeremy L Thompson DM dm; 188ccaff030SJeremy L Thompson DM dmviz; 189ccaff030SJeremy L Thompson Mat interpviz; 190ccaff030SJeremy L Thompson Ceed ceed; 191ccaff030SJeremy L Thompson Units units; 192ccaff030SJeremy L Thompson CeedVector qceed, qdotceed, gceed; 193ccaff030SJeremy L Thompson CeedOperator op_rhs, op_ifunction; 194ccaff030SJeremy L Thompson Vec M; 195ccaff030SJeremy L Thompson char outputfolder[PETSC_MAX_PATH_LEN]; 196ccaff030SJeremy L Thompson PetscInt contsteps; 197ccaff030SJeremy L Thompson }; 198ccaff030SJeremy L Thompson 199ccaff030SJeremy L Thompson struct Units_ { 200ccaff030SJeremy L Thompson // fundamental units 201ccaff030SJeremy L Thompson PetscScalar meter; 202ccaff030SJeremy L Thompson PetscScalar kilogram; 203ccaff030SJeremy L Thompson PetscScalar second; 204ccaff030SJeremy L Thompson PetscScalar Kelvin; 205ccaff030SJeremy L Thompson // derived units 206ccaff030SJeremy L Thompson PetscScalar Pascal; 207ccaff030SJeremy L Thompson PetscScalar JperkgK; 208ccaff030SJeremy L Thompson PetscScalar mpersquareds; 209ccaff030SJeremy L Thompson PetscScalar WpermK; 210ccaff030SJeremy L Thompson PetscScalar kgpercubicm; 211ccaff030SJeremy L Thompson PetscScalar kgpersquaredms; 212ccaff030SJeremy L Thompson PetscScalar Joulepercubicm; 213ccaff030SJeremy L Thompson }; 214ccaff030SJeremy L Thompson 215ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC; 216ccaff030SJeremy L Thompson struct SimpleBC_ { 217ccaff030SJeremy L Thompson PetscInt nwall, nslip[3]; 218c063f476Svaleriabarra PetscInt walls[6], slips[3][6]; 21907af6069Svaleriabarra PetscBool userbc; 220ccaff030SJeremy L Thompson }; 221ccaff030SJeremy L Thompson 222ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1). 223ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) { 224ccaff030SJeremy L Thompson return i >= 0 ? i : -(i+1); 225ccaff030SJeremy L Thompson } 226ccaff030SJeremy L Thompson 227ccaff030SJeremy L Thompson // Utility function to create local CEED restriction 228ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 229ccaff030SJeremy L Thompson CeedElemRestriction *Erestrict) { 230ccaff030SJeremy L Thompson 231ccaff030SJeremy L Thompson PetscSection section; 232ccaff030SJeremy L Thompson PetscInt c, cStart, cEnd, Nelem, Ndof, *erestrict, eoffset, nfields, dim; 233ccaff030SJeremy L Thompson PetscErrorCode ierr; 234ccaff030SJeremy L Thompson Vec Uloc; 235ccaff030SJeremy L Thompson 236ccaff030SJeremy L Thompson PetscFunctionBeginUser; 237ccaff030SJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 238ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm,§ion); CHKERRQ(ierr); 239ccaff030SJeremy L Thompson ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 240ccaff030SJeremy L Thompson PetscInt ncomp[nfields], fieldoff[nfields+1]; 241ccaff030SJeremy L Thompson fieldoff[0] = 0; 242ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 243ccaff030SJeremy L Thompson ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 244ccaff030SJeremy L Thompson fieldoff[f+1] = fieldoff[f] + ncomp[f]; 245ccaff030SJeremy L Thompson } 246ccaff030SJeremy L Thompson 247ccaff030SJeremy L Thompson ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 248ccaff030SJeremy L Thompson Nelem = cEnd - cStart; 249ccaff030SJeremy L Thompson ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 250ccaff030SJeremy L Thompson for (c=cStart,eoffset=0; c<cEnd; c++) { 251ccaff030SJeremy L Thompson PetscInt numindices, *indices, nnodes; 252ccaff030SJeremy L Thompson ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices, 253ccaff030SJeremy L Thompson &indices, NULL); CHKERRQ(ierr); 254ccaff030SJeremy L Thompson if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 255ccaff030SJeremy L Thompson PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 256ccaff030SJeremy L Thompson c); 257ccaff030SJeremy L Thompson nnodes = numindices / fieldoff[nfields]; 258ccaff030SJeremy L Thompson for (PetscInt i=0; i<nnodes; i++) { 259ccaff030SJeremy L Thompson // Check that indices are blocked by node and thus can be coalesced as a single field with 260ccaff030SJeremy L Thompson // fieldoff[nfields] = sum(ncomp) components. 261ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 262ccaff030SJeremy L Thompson for (PetscInt j=0; j<ncomp[f]; j++) { 263ccaff030SJeremy L Thompson if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j]) 264ccaff030SJeremy L Thompson != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j) 265ccaff030SJeremy L Thompson SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 266ccaff030SJeremy L Thompson "Cell %D closure indices not interlaced for node %D field %D component %D", 267ccaff030SJeremy L Thompson c, i, f, j); 268ccaff030SJeremy L Thompson } 269ccaff030SJeremy L Thompson } 270ccaff030SJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 271ccaff030SJeremy L Thompson PetscInt loc = Involute(indices[i*ncomp[0]]); 272*6f55dfd5Svaleriabarra erestrict[eoffset++] = loc; 273ccaff030SJeremy L Thompson } 274ccaff030SJeremy L Thompson ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices, 275ccaff030SJeremy L Thompson &indices, NULL); CHKERRQ(ierr); 276ccaff030SJeremy L Thompson } 277ccaff030SJeremy L Thompson if (eoffset != Nelem*PetscPowInt(P, dim)) SETERRQ3(PETSC_COMM_SELF, 278ccaff030SJeremy L Thompson PETSC_ERR_LIB, "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 279ccaff030SJeremy L Thompson PetscPowInt(P, dim),eoffset); 280ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 281ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 282ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 283*6f55dfd5Svaleriabarra CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields], 284*6f55dfd5Svaleriabarra 1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 285*6f55dfd5Svaleriabarra Erestrict); 286ccaff030SJeremy L Thompson ierr = PetscFree(erestrict); CHKERRQ(ierr); 287ccaff030SJeremy L Thompson PetscFunctionReturn(0); 288ccaff030SJeremy L Thompson } 289ccaff030SJeremy L Thompson 290ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 291ccaff030SJeremy L Thompson PetscErrorCode ierr; 292ccaff030SJeremy L Thompson PetscInt m; 293ccaff030SJeremy L Thompson 294ccaff030SJeremy L Thompson PetscFunctionBeginUser; 295ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 296ccaff030SJeremy L Thompson ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 297ccaff030SJeremy L Thompson PetscFunctionReturn(0); 298ccaff030SJeremy L Thompson } 299ccaff030SJeremy L Thompson 300ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) { 301ccaff030SJeremy L Thompson PetscErrorCode ierr; 302ccaff030SJeremy L Thompson PetscInt mceed,mpetsc; 303ccaff030SJeremy L Thompson PetscScalar *a; 304ccaff030SJeremy L Thompson 305ccaff030SJeremy L Thompson PetscFunctionBeginUser; 306ccaff030SJeremy L Thompson ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 307ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 308ccaff030SJeremy L Thompson if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 309380e0a58Svaleriabarra "Cannot place PETSc Vec of length %D in CeedVector of length %D", 310380e0a58Svaleriabarra mpetsc, mceed); 311ccaff030SJeremy L Thompson ierr = VecGetArray(p, &a); CHKERRQ(ierr); 312ccaff030SJeremy L Thompson CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 313ccaff030SJeremy L Thompson PetscFunctionReturn(0); 314ccaff030SJeremy L Thompson } 315ccaff030SJeremy L Thompson 316ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 317ccaff030SJeremy L Thompson PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 318ccaff030SJeremy L Thompson Vec cellGeomFVM, Vec gradFVM) { 319ccaff030SJeremy L Thompson PetscErrorCode ierr; 320ccaff030SJeremy L Thompson Vec Qbc; 321ccaff030SJeremy L Thompson 322ccaff030SJeremy L Thompson PetscFunctionBegin; 323ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 324ccaff030SJeremy L Thompson ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 325ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 326ccaff030SJeremy L Thompson PetscFunctionReturn(0); 327ccaff030SJeremy L Thompson } 328ccaff030SJeremy L Thompson 329ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u) 330ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G 331ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 332ccaff030SJeremy L Thompson PetscErrorCode ierr; 333ccaff030SJeremy L Thompson User user = *(User *)userData; 334ccaff030SJeremy L Thompson PetscScalar *q, *g; 335ccaff030SJeremy L Thompson Vec Qloc, Gloc; 336ccaff030SJeremy L Thompson 337ccaff030SJeremy L Thompson // Global-to-local 338ccaff030SJeremy L Thompson PetscFunctionBeginUser; 339ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 340ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 341ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 342ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 343ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 344ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 345ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 346ccaff030SJeremy L Thompson 347ccaff030SJeremy L Thompson // Ceed Vectors 348ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 349ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 350ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 351ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 352ccaff030SJeremy L Thompson 353ccaff030SJeremy L Thompson // Apply CEED operator 354ccaff030SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 355ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 356ccaff030SJeremy L Thompson 357ccaff030SJeremy L Thompson // Restore vectors 358ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 359ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 360ccaff030SJeremy L Thompson 361ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 362ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 363ccaff030SJeremy L Thompson 364ccaff030SJeremy L Thompson // Inverse of the lumped mass matrix 365ccaff030SJeremy L Thompson ierr = VecPointwiseMult(G, G, user->M); // M is Minv 366ccaff030SJeremy L Thompson CHKERRQ(ierr); 367ccaff030SJeremy L Thompson 368ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 369ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 370ccaff030SJeremy L Thompson PetscFunctionReturn(0); 371ccaff030SJeremy L Thompson } 372ccaff030SJeremy L Thompson 373ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 374ccaff030SJeremy L Thompson void *userData) { 375ccaff030SJeremy L Thompson PetscErrorCode ierr; 376ccaff030SJeremy L Thompson User user = *(User *)userData; 377ccaff030SJeremy L Thompson const PetscScalar *q, *qdot; 378ccaff030SJeremy L Thompson PetscScalar *g; 379ccaff030SJeremy L Thompson Vec Qloc, Qdotloc, Gloc; 380ccaff030SJeremy L Thompson 381ccaff030SJeremy L Thompson // Global-to-local 382ccaff030SJeremy L Thompson PetscFunctionBeginUser; 383ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 384ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 385ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 386ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 387ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 388ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 389ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 390ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 391ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 392ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 393ccaff030SJeremy L Thompson 394ccaff030SJeremy L Thompson // Ceed Vectors 395ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 396ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 397ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 398ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 399ccaff030SJeremy L Thompson (PetscScalar *)q); 400ccaff030SJeremy L Thompson CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 401ccaff030SJeremy L Thompson (PetscScalar *)qdot); 402ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 403ccaff030SJeremy L Thompson 404ccaff030SJeremy L Thompson // Apply CEED operator 405ccaff030SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 406ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 407ccaff030SJeremy L Thompson 408ccaff030SJeremy L Thompson // Restore vectors 409ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 410ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 411ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 412ccaff030SJeremy L Thompson 413ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 414ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 415ccaff030SJeremy L Thompson 416ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 417ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 418ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 419ccaff030SJeremy L Thompson PetscFunctionReturn(0); 420ccaff030SJeremy L Thompson } 421ccaff030SJeremy L Thompson 422ccaff030SJeremy L Thompson // User provided TS Monitor 423ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 424ccaff030SJeremy L Thompson Vec Q, void *ctx) { 425ccaff030SJeremy L Thompson User user = ctx; 426ccaff030SJeremy L Thompson Vec Qloc; 427ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 428ccaff030SJeremy L Thompson PetscViewer viewer; 429ccaff030SJeremy L Thompson PetscErrorCode ierr; 430ccaff030SJeremy L Thompson 431ccaff030SJeremy L Thompson // Set up output 432ccaff030SJeremy L Thompson PetscFunctionBeginUser; 433ccaff030SJeremy L Thompson // Print every 'outputfreq' steps 434ccaff030SJeremy L Thompson if (stepno % user->outputfreq != 0) 435ccaff030SJeremy L Thompson PetscFunctionReturn(0); 436ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 437ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 438ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 439ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 440ccaff030SJeremy L Thompson 441ccaff030SJeremy L Thompson // Output 442ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 443ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 444ccaff030SJeremy L Thompson CHKERRQ(ierr); 445ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 446ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 447ccaff030SJeremy L Thompson ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 448ccaff030SJeremy L Thompson if (user->dmviz) { 449ccaff030SJeremy L Thompson Vec Qrefined, Qrefined_loc; 450ccaff030SJeremy L Thompson char filepath_refined[PETSC_MAX_PATH_LEN]; 451ccaff030SJeremy L Thompson PetscViewer viewer_refined; 452ccaff030SJeremy L Thompson 453ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 454ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 455ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 456ccaff030SJeremy L Thompson CHKERRQ(ierr); 457ccaff030SJeremy L Thompson ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 458ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 459ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 460ccaff030SJeremy L Thompson CHKERRQ(ierr); 461ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 462ccaff030SJeremy L Thompson "%s/nsrefined-%03D.vtu", 463ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 464ccaff030SJeremy L Thompson CHKERRQ(ierr); 465ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 466ccaff030SJeremy L Thompson filepath_refined, 467ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 468ccaff030SJeremy L Thompson ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 469ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 470ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 471ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 472ccaff030SJeremy L Thompson } 473ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 474ccaff030SJeremy L Thompson 475ccaff030SJeremy L Thompson // Save data in a binary file for continuation of simulations 476ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 477ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 478ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 479ccaff030SJeremy L Thompson CHKERRQ(ierr); 480ccaff030SJeremy L Thompson ierr = VecView(Q, viewer); CHKERRQ(ierr); 481ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 482ccaff030SJeremy L Thompson 483ccaff030SJeremy L Thompson // Save time stamp 484ccaff030SJeremy L Thompson // Dimensionalize time back 485ccaff030SJeremy L Thompson time /= user->units->second; 486ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 487ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 488ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 489ccaff030SJeremy L Thompson CHKERRQ(ierr); 490ccaff030SJeremy L Thompson #if PETSC_VERSION_GE(3,13,0) 491ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 492ccaff030SJeremy L Thompson #else 493ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 494ccaff030SJeremy L Thompson #endif 495ccaff030SJeremy L Thompson CHKERRQ(ierr); 496ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 497ccaff030SJeremy L Thompson 498ccaff030SJeremy L Thompson PetscFunctionReturn(0); 499ccaff030SJeremy L Thompson } 500ccaff030SJeremy L Thompson 50127dd567eSvaleriabarra static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 502ccaff030SJeremy L Thompson CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 503ccaff030SJeremy L Thompson CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) { 504ccaff030SJeremy L Thompson PetscErrorCode ierr; 505ccaff030SJeremy L Thompson CeedVector multlvec; 506ccaff030SJeremy L Thompson Vec Multiplicity, MultiplicityLoc; 507ccaff030SJeremy L Thompson 508ccaff030SJeremy L Thompson ctxSetup->time = time; 509ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 510ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 511ccaff030SJeremy L Thompson CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 512ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 513ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 514ccaff030SJeremy L Thompson 515ccaff030SJeremy L Thompson // Fix multiplicity for output of ICs 516ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 517ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 518ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 519ccaff030SJeremy L Thompson CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 520ccaff030SJeremy L Thompson CeedVectorDestroy(&multlvec); 521ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 522ccaff030SJeremy L Thompson ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 523ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 524ccaff030SJeremy L Thompson CHKERRQ(ierr); 525ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 526ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 527ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 528ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 529ccaff030SJeremy L Thompson 530ccaff030SJeremy L Thompson PetscFunctionReturn(0); 531ccaff030SJeremy L Thompson } 532ccaff030SJeremy L Thompson 533ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 534ccaff030SJeremy L Thompson CeedElemRestriction restrictq, CeedBasis basisq, 535ccaff030SJeremy L Thompson CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 536ccaff030SJeremy L Thompson PetscErrorCode ierr; 537ccaff030SJeremy L Thompson CeedQFunction qf_mass; 538ccaff030SJeremy L Thompson CeedOperator op_mass; 539ccaff030SJeremy L Thompson CeedVector mceed; 540ccaff030SJeremy L Thompson Vec Mloc; 541ccaff030SJeremy L Thompson CeedInt ncompq, qdatasize; 542ccaff030SJeremy L Thompson 543ccaff030SJeremy L Thompson PetscFunctionBeginUser; 544ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 545ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 546ccaff030SJeremy L Thompson // Create the Q-function that defines the action of the mass operator 547ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 548ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 549ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 550ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 551ccaff030SJeremy L Thompson 552ccaff030SJeremy L Thompson // Create the mass operator 553ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 554ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 555ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", restrictqdi, 556ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 557ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 558ccaff030SJeremy L Thompson 559ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 560ccaff030SJeremy L Thompson ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 561ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 562ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 563ccaff030SJeremy L Thompson 564ccaff030SJeremy L Thompson { 565ccaff030SJeremy L Thompson // Compute a lumped mass matrix 566ccaff030SJeremy L Thompson CeedVector onesvec; 567ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 568ccaff030SJeremy L Thompson CeedVectorSetValue(onesvec, 1.0); 569ccaff030SJeremy L Thompson CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 570ccaff030SJeremy L Thompson CeedVectorDestroy(&onesvec); 571ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_mass); 572ccaff030SJeremy L Thompson CeedVectorDestroy(&mceed); 573ccaff030SJeremy L Thompson } 574ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_mass); 575ccaff030SJeremy L Thompson 576ccaff030SJeremy L Thompson ierr = VecZeroEntries(M); CHKERRQ(ierr); 577ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 578ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 579ccaff030SJeremy L Thompson 580ccaff030SJeremy L Thompson // Invert diagonally lumped mass vector for RHS function 581ccaff030SJeremy L Thompson ierr = VecReciprocal(M); CHKERRQ(ierr); 582ccaff030SJeremy L Thompson PetscFunctionReturn(0); 583ccaff030SJeremy L Thompson } 584ccaff030SJeremy L Thompson 585c063f476Svaleriabarra static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 586ff6701fcSJed Brown SimpleBC bc, void *ctxSetup) { 587ccaff030SJeremy L Thompson PetscErrorCode ierr; 588ccaff030SJeremy L Thompson 589ccaff030SJeremy L Thompson PetscFunctionBeginUser; 590ccaff030SJeremy L Thompson { 591ccaff030SJeremy L Thompson // Configure the finite element space and boundary conditions 592ccaff030SJeremy L Thompson PetscFE fe; 593ccaff030SJeremy L Thompson PetscInt ncompq = 5; 594ff6701fcSJed Brown ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 595ff6701fcSJed Brown PETSC_FALSE, degree, PETSC_DECIDE, 596ff6701fcSJed Brown &fe); 597ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 598ccaff030SJeremy L Thompson ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr); 599ccaff030SJeremy L Thompson ierr = DMCreateDS(dm); CHKERRQ(ierr); 60007af6069Svaleriabarra { 60107af6069Svaleriabarra PetscInt comps[1] = {1}; 60207af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 60307af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[0], 60407af6069Svaleriabarra bc->slips[0], ctxSetup); CHKERRQ(ierr); 60507af6069Svaleriabarra comps[0] = 2; 60607af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 60707af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[1], 60807af6069Svaleriabarra bc->slips[1], ctxSetup); CHKERRQ(ierr); 60907af6069Svaleriabarra comps[0] = 3; 61007af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 61107af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[2], 61207af6069Svaleriabarra bc->slips[2], ctxSetup); CHKERRQ(ierr); 61307af6069Svaleriabarra } 61407af6069Svaleriabarra if (bc->userbc == PETSC_TRUE) { 61507af6069Svaleriabarra for (PetscInt c = 0; c < 3; c++) { 61607af6069Svaleriabarra for (PetscInt s = 0; s < bc->nslip[c]; s++) { 61707af6069Svaleriabarra for (PetscInt w = 0; w < bc->nwall; w++) { 61807af6069Svaleriabarra if (bc->slips[c][s] == bc->walls[w]) 61907af6069Svaleriabarra SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 62007af6069Svaleriabarra "Boundary condition already set on face %D!\n", bc->walls[w]); 62107af6069Svaleriabarra 62207af6069Svaleriabarra } 62307af6069Svaleriabarra } 62407af6069Svaleriabarra } 62507af6069Svaleriabarra } 62635d2ed1bSvaleriabarra // Wall boundary conditions are zero energy density and zero flux for 62735d2ed1bSvaleriabarra // velocity in advection/advection2d, and zero velocity and zero flux 62835d2ed1bSvaleriabarra // for mass density and energy density in density_current 629ccaff030SJeremy L Thompson { 63035d2ed1bSvaleriabarra if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 63135d2ed1bSvaleriabarra PetscInt comps[1] = {4}; 63235d2ed1bSvaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 63335d2ed1bSvaleriabarra 1, comps, (void(*)(void))problem->bc, 63435d2ed1bSvaleriabarra bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 63535d2ed1bSvaleriabarra } else if (problem->bc == Exact_DC) { 636ccaff030SJeremy L Thompson PetscInt comps[3] = {1, 2, 3}; 637ccaff030SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 638ccaff030SJeremy L Thompson 3, comps, (void(*)(void))problem->bc, 639ccaff030SJeremy L Thompson bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 64007af6069Svaleriabarra } else 64107af6069Svaleriabarra SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, 642c063f476Svaleriabarra "Undefined boundary conditions for this problem"); 643ccaff030SJeremy L Thompson } 644ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 645ccaff030SJeremy L Thompson CHKERRQ(ierr); 646ccaff030SJeremy L Thompson ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 647ccaff030SJeremy L Thompson } 648ccaff030SJeremy L Thompson { 649ccaff030SJeremy L Thompson // Empty name for conserved field (because there is only one field) 650ccaff030SJeremy L Thompson PetscSection section; 651ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 652ccaff030SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 653ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 654ccaff030SJeremy L Thompson CHKERRQ(ierr); 655ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 656ccaff030SJeremy L Thompson CHKERRQ(ierr); 657ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 658ccaff030SJeremy L Thompson CHKERRQ(ierr); 659ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 660ccaff030SJeremy L Thompson CHKERRQ(ierr); 661ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 662ccaff030SJeremy L Thompson CHKERRQ(ierr); 663ccaff030SJeremy L Thompson } 664ccaff030SJeremy L Thompson PetscFunctionReturn(0); 665ccaff030SJeremy L Thompson } 666ccaff030SJeremy L Thompson 667ccaff030SJeremy L Thompson int main(int argc, char **argv) { 668ccaff030SJeremy L Thompson PetscInt ierr; 669ccaff030SJeremy L Thompson MPI_Comm comm; 670380e0a58Svaleriabarra DM dm, dmcoord, dmviz; 671ccaff030SJeremy L Thompson Mat interpviz; 672ccaff030SJeremy L Thompson TS ts; 673ccaff030SJeremy L Thompson TSAdapt adapt; 674ccaff030SJeremy L Thompson User user; 675ccaff030SJeremy L Thompson Units units; 676ccaff030SJeremy L Thompson char ceedresource[4096] = "/cpu/self"; 6770261bf83Svaleriabarra PetscInt cStart, cEnd, localNelem, lnodes, gnodes, steps; 678ccaff030SJeremy L Thompson const PetscInt ncompq = 5; 679ccaff030SJeremy L Thompson PetscMPIInt rank; 680ccaff030SJeremy L Thompson PetscScalar ftime; 681ccaff030SJeremy L Thompson Vec Q, Qloc, Xloc; 682ccaff030SJeremy L Thompson Ceed ceed; 683ccaff030SJeremy L Thompson CeedInt numP, numQ; 684ccaff030SJeremy L Thompson CeedVector xcorners, qdata, q0ceed; 685ccaff030SJeremy L Thompson CeedBasis basisx, basisxc, basisq; 686*6f55dfd5Svaleriabarra CeedElemRestriction restrictx, restrictq, restrictqdi; 687ccaff030SJeremy L Thompson CeedQFunction qf_setup, qf_ics, qf_rhs, qf_ifunction; 688ccaff030SJeremy L Thompson CeedOperator op_setup, op_ics; 689ccaff030SJeremy L Thompson CeedScalar Rd; 6900261bf83Svaleriabarra CeedMemType memtyperequested; 691ccaff030SJeremy L Thompson PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 692ccaff030SJeremy L Thompson kgpersquaredms, Joulepercubicm; 693ccaff030SJeremy L Thompson problemType problemChoice; 694ccaff030SJeremy L Thompson problemData *problem = NULL; 695ccaff030SJeremy L Thompson StabilizationType stab; 6968ca749c7Svaleriabarra testType testChoice; 6978ca749c7Svaleriabarra testData *test = NULL; 6988ca749c7Svaleriabarra PetscBool implicit; 699cb3e2689Svaleriabarra PetscInt viz_refine = 0; 700ccaff030SJeremy L Thompson struct SimpleBC_ bc = { 70107af6069Svaleriabarra .nslip = {2, 2, 2}, 70207af6069Svaleriabarra .slips = {{5, 6}, {3, 4}, {1, 2}} 703ccaff030SJeremy L Thompson }; 704ccaff030SJeremy L Thompson double start, cpu_time_used; 7050261bf83Svaleriabarra // Check PETSc CUDA support 7060261bf83Svaleriabarra PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 7070261bf83Svaleriabarra // *INDENT-OFF* 7080261bf83Svaleriabarra #ifdef PETSC_HAVE_CUDA 7090261bf83Svaleriabarra petschavecuda = PETSC_TRUE; 7100261bf83Svaleriabarra #else 7110261bf83Svaleriabarra petschavecuda = PETSC_FALSE; 7120261bf83Svaleriabarra #endif 7130261bf83Svaleriabarra // *INDENT-ON* 714ccaff030SJeremy L Thompson 715ccaff030SJeremy L Thompson // Create the libCEED contexts 716ccaff030SJeremy L Thompson PetscScalar meter = 1e-2; // 1 meter in scaled length units 717ccaff030SJeremy L Thompson PetscScalar second = 1e-2; // 1 second in scaled time units 718ccaff030SJeremy L Thompson PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 719ccaff030SJeremy L Thompson PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 720ccaff030SJeremy L Thompson CeedScalar theta0 = 300.; // K 721ccaff030SJeremy L Thompson CeedScalar thetaC = -15.; // K 722ccaff030SJeremy L Thompson CeedScalar P0 = 1.e5; // Pa 723ccaff030SJeremy L Thompson CeedScalar N = 0.01; // 1/s 724ccaff030SJeremy L Thompson CeedScalar cv = 717.; // J/(kg K) 725ccaff030SJeremy L Thompson CeedScalar cp = 1004.; // J/(kg K) 726ccaff030SJeremy L Thompson CeedScalar g = 9.81; // m/s^2 727ccaff030SJeremy L Thompson CeedScalar lambda = -2./3.; // - 728ccaff030SJeremy L Thompson CeedScalar mu = 75.; // Pa s, dynamic viscosity 729ccaff030SJeremy L Thompson // mu = 75 is not physical for air, but is good for numerical stability 730ccaff030SJeremy L Thompson CeedScalar k = 0.02638; // W/(m K) 731ccaff030SJeremy L Thompson CeedScalar CtauS = 0.; // dimensionless 732ccaff030SJeremy L Thompson CeedScalar strong_form = 0.; // [0,1] 733ccaff030SJeremy L Thompson PetscScalar lx = 8000.; // m 734ccaff030SJeremy L Thompson PetscScalar ly = 8000.; // m 735ccaff030SJeremy L Thompson PetscScalar lz = 4000.; // m 736ccaff030SJeremy L Thompson CeedScalar rc = 1000.; // m (Radius of bubble) 737ccaff030SJeremy L Thompson PetscScalar resx = 1000.; // m (resolution in x) 738ccaff030SJeremy L Thompson PetscScalar resy = 1000.; // m (resolution in y) 739ccaff030SJeremy L Thompson PetscScalar resz = 1000.; // m (resolution in z) 740ccaff030SJeremy L Thompson PetscInt outputfreq = 10; // - 741ccaff030SJeremy L Thompson PetscInt contsteps = 0; // - 742ff6701fcSJed Brown PetscInt degree = 1; // - 743ccaff030SJeremy L Thompson PetscInt qextra = 2; // - 744ccaff030SJeremy L Thompson PetscReal center[3], dc_axis[3] = {0, 0, 0}; 745ccaff030SJeremy L Thompson 746ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 747ccaff030SJeremy L Thompson if (ierr) return ierr; 748ccaff030SJeremy L Thompson 749ccaff030SJeremy L Thompson // Allocate PETSc context 750ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 751ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 752ccaff030SJeremy L Thompson 753ccaff030SJeremy L Thompson // Parse command line options 754ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 755ccaff030SJeremy L Thompson ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 756ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 757ccaff030SJeremy L Thompson ierr = PetscOptionsString("-ceed", "CEED resource specifier", 758ccaff030SJeremy L Thompson NULL, ceedresource, ceedresource, 759ccaff030SJeremy L Thompson sizeof(ceedresource), NULL); CHKERRQ(ierr); 7608ca749c7Svaleriabarra testChoice = TEST_NONE; 7618ca749c7Svaleriabarra ierr = PetscOptionsEnum("-test", "Run tests", NULL, 7628ca749c7Svaleriabarra testTypes, (PetscEnum)testChoice, 7638ca749c7Svaleriabarra (PetscEnum *)&testChoice, 7648ca749c7Svaleriabarra NULL); CHKERRQ(ierr); 7658ca749c7Svaleriabarra test = &testOptions[testChoice]; 766ccaff030SJeremy L Thompson problemChoice = NS_DENSITY_CURRENT; 767ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 768ccaff030SJeremy L Thompson problemTypes, (PetscEnum)problemChoice, 769ccaff030SJeremy L Thompson (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 770ccaff030SJeremy L Thompson problem = &problemOptions[problemChoice]; 771ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 772ccaff030SJeremy L Thompson StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 773ccaff030SJeremy L Thompson (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 774ccaff030SJeremy L Thompson ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 775ccaff030SJeremy L Thompson NULL, implicit=PETSC_FALSE, &implicit, NULL); 776ccaff030SJeremy L Thompson CHKERRQ(ierr); 7778c662231Svaleriabarra if (!implicit && stab != STAB_NONE) { 7788c662231Svaleriabarra ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 7798c662231Svaleriabarra CHKERRQ(ierr); 7808c662231Svaleriabarra } 781ccaff030SJeremy L Thompson { 782ccaff030SJeremy L Thompson PetscInt len; 783ccaff030SJeremy L Thompson PetscBool flg; 784ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray("-bc_wall", 785ccaff030SJeremy L Thompson "Use wall boundary conditions on this list of faces", 786ccaff030SJeremy L Thompson NULL, bc.walls, 787ccaff030SJeremy L Thompson (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 788ccaff030SJeremy L Thompson &len), &flg); CHKERRQ(ierr); 789ccaff030SJeremy L Thompson if (flg) bc.nwall = len; 790ccaff030SJeremy L Thompson for (PetscInt j=0; j<3; j++) { 791ccaff030SJeremy L Thompson const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 792ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray(flags[j], 793ccaff030SJeremy L Thompson "Use slip boundary conditions on this list of faces", 794ccaff030SJeremy L Thompson NULL, bc.slips[j], 795ccaff030SJeremy L Thompson (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 796ccaff030SJeremy L Thompson &len), &flg); 797ccaff030SJeremy L Thompson CHKERRQ(ierr); 79807af6069Svaleriabarra if (flg) { 79907af6069Svaleriabarra bc.nslip[j] = len; 80007af6069Svaleriabarra bc.userbc = PETSC_TRUE; 80107af6069Svaleriabarra } 802ccaff030SJeremy L Thompson } 803ccaff030SJeremy L Thompson } 804cb3e2689Svaleriabarra ierr = PetscOptionsInt("-viz_refine", 805cb3e2689Svaleriabarra "Regular refinement levels for visualization", 806cb3e2689Svaleriabarra NULL, viz_refine, &viz_refine, NULL); 807ccaff030SJeremy L Thompson CHKERRQ(ierr); 808ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 809ccaff030SJeremy L Thompson NULL, meter, &meter, NULL); CHKERRQ(ierr); 810ccaff030SJeremy L Thompson meter = fabs(meter); 811ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 812ccaff030SJeremy L Thompson NULL, second, &second, NULL); CHKERRQ(ierr); 813ccaff030SJeremy L Thompson second = fabs(second); 814ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 815ccaff030SJeremy L Thompson NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson kilogram = fabs(kilogram); 817ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_Kelvin", 818ccaff030SJeremy L Thompson "1 Kelvin in scaled temperature units", 819ccaff030SJeremy L Thompson NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson Kelvin = fabs(Kelvin); 821ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 822ccaff030SJeremy L Thompson NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 823ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 824ccaff030SJeremy L Thompson NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 825ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 826ccaff030SJeremy L Thompson NULL, P0, &P0, NULL); CHKERRQ(ierr); 827ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 828ccaff030SJeremy L Thompson NULL, N, &N, NULL); CHKERRQ(ierr); 829ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 830ccaff030SJeremy L Thompson NULL, cv, &cv, NULL); CHKERRQ(ierr); 831ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 832ccaff030SJeremy L Thompson NULL, cp, &cp, NULL); CHKERRQ(ierr); 833ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 834ccaff030SJeremy L Thompson NULL, g, &g, NULL); CHKERRQ(ierr); 835ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lambda", 836ccaff030SJeremy L Thompson "Stokes hypothesis second viscosity coefficient", 837ccaff030SJeremy L Thompson NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 838ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 839ccaff030SJeremy L Thompson NULL, mu, &mu, NULL); CHKERRQ(ierr); 840ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-k", "Thermal conductivity", 841ccaff030SJeremy L Thompson NULL, k, &k, NULL); CHKERRQ(ierr); 842ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-CtauS", 843ccaff030SJeremy L Thompson "Scale coefficient for tau (nondimensional)", 844ccaff030SJeremy L Thompson NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 84540c555b9Svaleriabarra if (stab == STAB_NONE && CtauS != 0) { 846c063f476Svaleriabarra ierr = PetscPrintf(comm, 847c063f476Svaleriabarra "Warning! Use -CtauS only with -stab su or -stab supg\n"); 84840c555b9Svaleriabarra CHKERRQ(ierr); 84940c555b9Svaleriabarra } 850ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-strong_form", 851ccaff030SJeremy L Thompson "Strong (1) or weak/integrated by parts (0) advection residual", 852ccaff030SJeremy L Thompson NULL, strong_form, &strong_form, NULL); 853ccaff030SJeremy L Thompson CHKERRQ(ierr); 854559ec21dSvaleriabarra if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) { 855559ec21dSvaleriabarra ierr = PetscPrintf(comm, 856559ec21dSvaleriabarra "Warning! Problem density_current does not support -CtauS or -strong_form\n"); 857559ec21dSvaleriabarra CHKERRQ(ierr); 858559ec21dSvaleriabarra } 859ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 860ccaff030SJeremy L Thompson NULL, lx, &lx, NULL); CHKERRQ(ierr); 861ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 862ccaff030SJeremy L Thompson NULL, ly, &ly, NULL); CHKERRQ(ierr); 863ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 864ccaff030SJeremy L Thompson NULL, lz, &lz, NULL); CHKERRQ(ierr); 865ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 866ccaff030SJeremy L Thompson NULL, rc, &rc, NULL); CHKERRQ(ierr); 867ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resx","Target resolution in x", 868ccaff030SJeremy L Thompson NULL, resx, &resx, NULL); CHKERRQ(ierr); 869ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resy","Target resolution in y", 870ccaff030SJeremy L Thompson NULL, resy, &resy, NULL); CHKERRQ(ierr); 871ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resz","Target resolution in z", 872ccaff030SJeremy L Thompson NULL, resz, &resz, NULL); CHKERRQ(ierr); 873ccaff030SJeremy L Thompson PetscInt n = problem->dim; 874ccaff030SJeremy L Thompson center[0] = 0.5 * lx; 875ccaff030SJeremy L Thompson center[1] = 0.5 * ly; 876ccaff030SJeremy L Thompson center[2] = 0.5 * lz; 877ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-center", "Location of bubble center", 878ccaff030SJeremy L Thompson NULL, center, &n, NULL); CHKERRQ(ierr); 879ccaff030SJeremy L Thompson n = problem->dim; 880ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-dc_axis", 881ccaff030SJeremy L Thompson "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 882ccaff030SJeremy L Thompson NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 883ccaff030SJeremy L Thompson { 884ccaff030SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 885ccaff030SJeremy L Thompson PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 886ccaff030SJeremy L Thompson if (norm > 0) { 887ccaff030SJeremy L Thompson for (int i=0; i<3; i++) dc_axis[i] /= norm; 888ccaff030SJeremy L Thompson } 889ccaff030SJeremy L Thompson } 890ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-output_freq", 891ccaff030SJeremy L Thompson "Frequency of output, in number of steps", 892ccaff030SJeremy L Thompson NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 893ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-continue", "Continue from previous solution", 894ccaff030SJeremy L Thompson NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 895ff6701fcSJed Brown ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 896ff6701fcSJed Brown NULL, degree, °ree, NULL); CHKERRQ(ierr); 897ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 898ccaff030SJeremy L Thompson NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 8993e4faa5aSvaleriabarra ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr); 900ccaff030SJeremy L Thompson ierr = PetscOptionsString("-of", "Output folder", 901ccaff030SJeremy L Thompson NULL, user->outputfolder, user->outputfolder, 902ccaff030SJeremy L Thompson sizeof(user->outputfolder), NULL); CHKERRQ(ierr); 9030261bf83Svaleriabarra memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 9040261bf83Svaleriabarra ierr = PetscOptionsEnum("-memtype", 9050261bf83Svaleriabarra "CEED MemType requested", NULL, 9060261bf83Svaleriabarra memTypes, (PetscEnum)memtyperequested, 9070261bf83Svaleriabarra (PetscEnum *)&memtyperequested, &setmemtyperequest); 9080261bf83Svaleriabarra CHKERRQ(ierr); 909ccaff030SJeremy L Thompson ierr = PetscOptionsEnd(); CHKERRQ(ierr); 910ccaff030SJeremy L Thompson 911ccaff030SJeremy L Thompson // Define derived units 912ccaff030SJeremy L Thompson Pascal = kilogram / (meter * PetscSqr(second)); 913ccaff030SJeremy L Thompson JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 914ccaff030SJeremy L Thompson mpersquareds = meter / PetscSqr(second); 915ccaff030SJeremy L Thompson WpermK = kilogram * meter / (pow(second,3) * Kelvin); 916ccaff030SJeremy L Thompson kgpercubicm = kilogram / pow(meter,3); 917ccaff030SJeremy L Thompson kgpersquaredms = kilogram / (PetscSqr(meter) * second); 918ccaff030SJeremy L Thompson Joulepercubicm = kilogram / (meter * PetscSqr(second)); 919ccaff030SJeremy L Thompson 920ccaff030SJeremy L Thompson // Scale variables to desired units 921ccaff030SJeremy L Thompson theta0 *= Kelvin; 922ccaff030SJeremy L Thompson thetaC *= Kelvin; 923ccaff030SJeremy L Thompson P0 *= Pascal; 924ccaff030SJeremy L Thompson N *= (1./second); 925ccaff030SJeremy L Thompson cv *= JperkgK; 926ccaff030SJeremy L Thompson cp *= JperkgK; 927ccaff030SJeremy L Thompson Rd = cp - cv; 928ccaff030SJeremy L Thompson g *= mpersquareds; 929ccaff030SJeremy L Thompson mu *= Pascal * second; 930ccaff030SJeremy L Thompson k *= WpermK; 931ccaff030SJeremy L Thompson lx = fabs(lx) * meter; 932ccaff030SJeremy L Thompson ly = fabs(ly) * meter; 933ccaff030SJeremy L Thompson lz = fabs(lz) * meter; 934ccaff030SJeremy L Thompson rc = fabs(rc) * meter; 935ccaff030SJeremy L Thompson resx = fabs(resx) * meter; 936ccaff030SJeremy L Thompson resy = fabs(resy) * meter; 937ccaff030SJeremy L Thompson resz = fabs(resz) * meter; 938ccaff030SJeremy L Thompson for (int i=0; i<3; i++) center[i] *= meter; 939ccaff030SJeremy L Thompson 940ccaff030SJeremy L Thompson const CeedInt dim = problem->dim, ncompx = problem->dim, 941ccaff030SJeremy L Thompson qdatasize = problem->qdatasize; 942ccaff030SJeremy L Thompson // Set up the libCEED context 943ccaff030SJeremy L Thompson struct SetupContext_ ctxSetup = { 944ccaff030SJeremy L Thompson .theta0 = theta0, 945ccaff030SJeremy L Thompson .thetaC = thetaC, 946ccaff030SJeremy L Thompson .P0 = P0, 947ccaff030SJeremy L Thompson .N = N, 948ccaff030SJeremy L Thompson .cv = cv, 949ccaff030SJeremy L Thompson .cp = cp, 950ccaff030SJeremy L Thompson .Rd = Rd, 951ccaff030SJeremy L Thompson .g = g, 952ccaff030SJeremy L Thompson .rc = rc, 953ccaff030SJeremy L Thompson .lx = lx, 954ccaff030SJeremy L Thompson .ly = ly, 955ccaff030SJeremy L Thompson .lz = lz, 956ccaff030SJeremy L Thompson .center[0] = center[0], 957ccaff030SJeremy L Thompson .center[1] = center[1], 958ccaff030SJeremy L Thompson .center[2] = center[2], 959ccaff030SJeremy L Thompson .dc_axis[0] = dc_axis[0], 960ccaff030SJeremy L Thompson .dc_axis[1] = dc_axis[1], 961ccaff030SJeremy L Thompson .dc_axis[2] = dc_axis[2], 962ccaff030SJeremy L Thompson .time = 0, 963ccaff030SJeremy L Thompson }; 964ccaff030SJeremy L Thompson 965380e0a58Svaleriabarra // Create the mesh 966ccaff030SJeremy L Thompson { 967ccaff030SJeremy L Thompson const PetscReal scale[3] = {lx, ly, lz}; 968ccaff030SJeremy L Thompson ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 969ff62f740Svaleriabarra NULL, PETSC_TRUE, &dm); 970ccaff030SJeremy L Thompson CHKERRQ(ierr); 971ccaff030SJeremy L Thompson } 972380e0a58Svaleriabarra 973380e0a58Svaleriabarra // Distribute the mesh over processes 974380e0a58Svaleriabarra { 975ccaff030SJeremy L Thompson DM dmDist = NULL; 976ccaff030SJeremy L Thompson PetscPartitioner part; 977ccaff030SJeremy L Thompson 978ccaff030SJeremy L Thompson ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 979ccaff030SJeremy L Thompson ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 980ccaff030SJeremy L Thompson ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 981ccaff030SJeremy L Thompson if (dmDist) { 982ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 983ccaff030SJeremy L Thompson dm = dmDist; 984ccaff030SJeremy L Thompson } 985ccaff030SJeremy L Thompson } 986ccaff030SJeremy L Thompson ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 987ccaff030SJeremy L Thompson 988380e0a58Svaleriabarra // Setup DM 989ccaff030SJeremy L Thompson ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 990ccaff030SJeremy L Thompson ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 991ff6701fcSJed Brown ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr); 992380e0a58Svaleriabarra 993380e0a58Svaleriabarra // Refine DM for high-order viz 994ccaff030SJeremy L Thompson dmviz = NULL; 995ccaff030SJeremy L Thompson interpviz = NULL; 996ccaff030SJeremy L Thompson if (viz_refine) { 997ff6701fcSJed Brown DM dmhierarchy[viz_refine+1]; 998ff6701fcSJed Brown 999ccaff030SJeremy L Thompson ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 1000ff6701fcSJed Brown dmhierarchy[0] = dm; 1001ff6701fcSJed Brown for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 1002ff6701fcSJed Brown Mat interp_next; 1003ff6701fcSJed Brown 1004ff6701fcSJed Brown ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 1005ccaff030SJeremy L Thompson CHKERRQ(ierr); 1006ff6701fcSJed Brown ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 1007ff6701fcSJed Brown d = (d + 1) / 2; 1008ff6701fcSJed Brown if (i + 1 == viz_refine) d = 1; 1009ff6701fcSJed Brown ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr); 1010ff6701fcSJed Brown ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 1011ff6701fcSJed Brown &interp_next, NULL); CHKERRQ(ierr); 1012ff6701fcSJed Brown if (!i) interpviz = interp_next; 1013ff6701fcSJed Brown else { 1014ff6701fcSJed Brown Mat C; 1015ff6701fcSJed Brown ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 1016ff6701fcSJed Brown PETSC_DECIDE, &C); CHKERRQ(ierr); 1017ff6701fcSJed Brown ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 1018ff6701fcSJed Brown ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1019ff6701fcSJed Brown interpviz = C; 1020ff6701fcSJed Brown } 1021ff6701fcSJed Brown } 1022cb3e2689Svaleriabarra for (PetscInt i=1; i<viz_refine; i++) { 1023ff6701fcSJed Brown ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 1024cb3e2689Svaleriabarra } 1025ff6701fcSJed Brown dmviz = dmhierarchy[viz_refine]; 1026ccaff030SJeremy L Thompson } 1027ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 1028ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 1029ccaff030SJeremy L Thompson ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 1030ccaff030SJeremy L Thompson lnodes /= ncompq; 1031ccaff030SJeremy L Thompson 10320261bf83Svaleriabarra // Initialize CEED 10330261bf83Svaleriabarra CeedInit(ceedresource, &ceed); 10340261bf83Svaleriabarra // Set memtype 10350261bf83Svaleriabarra CeedMemType memtypebackend; 10360261bf83Svaleriabarra CeedGetPreferredMemType(ceed, &memtypebackend); 10370261bf83Svaleriabarra // Check memtype compatibility 10380261bf83Svaleriabarra if (!setmemtyperequest) 10390261bf83Svaleriabarra memtyperequested = memtypebackend; 10400261bf83Svaleriabarra else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 10410261bf83Svaleriabarra SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 10420261bf83Svaleriabarra "PETSc was not built with CUDA. " 10430261bf83Svaleriabarra "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 10440261bf83Svaleriabarra 10450261bf83Svaleriabarra // Set number of 1D nodes and quadrature points 10460261bf83Svaleriabarra numP = degree + 1; 10470261bf83Svaleriabarra numQ = numP + qextra; 10480261bf83Svaleriabarra 10490261bf83Svaleriabarra // Print summary 10500261bf83Svaleriabarra if (testChoice == TEST_NONE) { 1051ccaff030SJeremy L Thompson CeedInt gdofs, odofs; 1052ccaff030SJeremy L Thompson int comm_size; 1053ccaff030SJeremy L Thompson char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 1054ccaff030SJeremy L Thompson ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 1055ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 10560261bf83Svaleriabarra gnodes = gdofs/ncompq; 1057ccaff030SJeremy L Thompson ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 1058ccaff030SJeremy L Thompson ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 1059ccaff030SJeremy L Thompson sizeof(box_faces_str), NULL); CHKERRQ(ierr); 10600261bf83Svaleriabarra const char *usedresource; 10610261bf83Svaleriabarra CeedGetResource(ceed, &usedresource); 10620261bf83Svaleriabarra 10630261bf83Svaleriabarra ierr = PetscPrintf(comm, 1064*6f55dfd5Svaleriabarra "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 1065*6f55dfd5Svaleriabarra " rank(s) : %d\n" 1066*6f55dfd5Svaleriabarra " Problem to solve : %s\n" 1067*6f55dfd5Svaleriabarra " Stabilization : %s\n" 10680261bf83Svaleriabarra " PETSc:\n" 1069*6f55dfd5Svaleriabarra " Box Faces: : %s\n" 10700261bf83Svaleriabarra " libCEED:\n" 10710261bf83Svaleriabarra " libCEED Backend : %s\n" 10720261bf83Svaleriabarra " libCEED Backend MemType : %s\n" 10730261bf83Svaleriabarra " libCEED User Requested MemType : %s\n" 10740261bf83Svaleriabarra " Mesh:\n" 10750261bf83Svaleriabarra " Number of 1D Basis Nodes (P) : %d\n" 10760261bf83Svaleriabarra " Number of 1D Quadrature Points (Q) : %d\n" 10770261bf83Svaleriabarra " Global DoFs : %D\n" 10780261bf83Svaleriabarra " Owned DoFs : %D\n" 10790261bf83Svaleriabarra " DoFs per node : %D\n" 10800261bf83Svaleriabarra " Global nodes : %D\n" 10810261bf83Svaleriabarra " Owned nodes : %D\n", 1082*6f55dfd5Svaleriabarra comm_size, problemTypes[problemChoice], 1083*6f55dfd5Svaleriabarra StabilizationTypes[stab], box_faces_str, usedresource, 10840261bf83Svaleriabarra CeedMemTypes[memtypebackend], 10850261bf83Svaleriabarra (setmemtyperequest) ? 10860261bf83Svaleriabarra CeedMemTypes[memtyperequested] : "none", 10870261bf83Svaleriabarra numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes); 1088ccaff030SJeremy L Thompson CHKERRQ(ierr); 1089ccaff030SJeremy L Thompson } 1090ccaff030SJeremy L Thompson 1091ccaff030SJeremy L Thompson // Set up global mass vector 1092ccaff030SJeremy L Thompson ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 1093ccaff030SJeremy L Thompson 10940261bf83Svaleriabarra // Set up libCEED 1095ccaff030SJeremy L Thompson // CEED Bases 1096ccaff030SJeremy L Thompson CeedInit(ceedresource, &ceed); 1097ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 1098ccaff030SJeremy L Thompson &basisq); 1099ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 1100ccaff030SJeremy L Thompson &basisx); 1101ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 1102ccaff030SJeremy L Thompson CEED_GAUSS_LOBATTO, &basisxc); 1103ccaff030SJeremy L Thompson 1104ccaff030SJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 1105ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 1106ccaff030SJeremy L Thompson CHKERRQ(ierr); 1107ccaff030SJeremy L Thompson 1108ccaff030SJeremy L Thompson // CEED Restrictions 1109ccaff030SJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, degree+1, &restrictq); 1110ccaff030SJeremy L Thompson CHKERRQ(ierr); 1111ccaff030SJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, &restrictx); CHKERRQ(ierr); 1112ccaff030SJeremy L Thompson DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 1113ccaff030SJeremy L Thompson localNelem = cEnd - cStart; 1114ccaff030SJeremy L Thompson CeedInt numQdim = CeedIntPow(numQ, dim); 1115ccaff030SJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, localNelem, numQdim, 1116*6f55dfd5Svaleriabarra qdatasize, qdatasize*localNelem*numQdim, 1117ccaff030SJeremy L Thompson CEED_STRIDES_BACKEND, &restrictqdi); 1118ccaff030SJeremy L Thompson 1119ccaff030SJeremy L Thompson ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1120ccaff030SJeremy L Thompson ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1121ccaff030SJeremy L Thompson 1122ccaff030SJeremy L Thompson // Create the CEED vectors that will be needed in setup 1123ccaff030SJeremy L Thompson CeedInt Nqpts; 1124ccaff030SJeremy L Thompson CeedBasisGetNumQuadraturePoints(basisq, &Nqpts); 1125ccaff030SJeremy L Thompson CeedVectorCreate(ceed, qdatasize*localNelem*Nqpts, &qdata); 1126ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1127ccaff030SJeremy L Thompson 1128ccaff030SJeremy L Thompson // Create the Q-function that builds the quadrature data for the NS operator 1129ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->setup, problem->setup_loc, 1130ccaff030SJeremy L Thompson &qf_setup); 1131ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_setup, "dx", ncompx*dim, CEED_EVAL_GRAD); 1132ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 1133ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_setup, "qdata", qdatasize, CEED_EVAL_NONE); 1134ccaff030SJeremy L Thompson 1135ccaff030SJeremy L Thompson // Create the Q-function that sets the ICs of the operator 1136ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1137ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1138ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1139ccaff030SJeremy L Thompson 1140ccaff030SJeremy L Thompson qf_rhs = NULL; 1141ccaff030SJeremy L Thompson if (problem->apply_rhs) { // Create the Q-function that defines the action of the RHS operator 1142ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->apply_rhs, 1143ccaff030SJeremy L Thompson problem->apply_rhs_loc, &qf_rhs); 1144ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "q", ncompq, CEED_EVAL_INTERP); 1145ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "dq", ncompq*dim, CEED_EVAL_GRAD); 1146ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "qdata", qdatasize, CEED_EVAL_NONE); 1147ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "x", ncompx, CEED_EVAL_INTERP); 1148ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_rhs, "v", ncompq, CEED_EVAL_INTERP); 1149ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_rhs, "dv", ncompq*dim, CEED_EVAL_GRAD); 1150ccaff030SJeremy L Thompson } 1151ccaff030SJeremy L Thompson 1152ccaff030SJeremy L Thompson qf_ifunction = NULL; 1153ccaff030SJeremy L Thompson if (problem->apply_ifunction) { // Create the Q-function that defines the action of the IFunction 1154ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->apply_ifunction, 1155ccaff030SJeremy L Thompson problem->apply_ifunction_loc, &qf_ifunction); 1156ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "q", ncompq, CEED_EVAL_INTERP); 1157ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "dq", ncompq*dim, CEED_EVAL_GRAD); 1158ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "qdot", ncompq, CEED_EVAL_INTERP); 1159ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "qdata", qdatasize, CEED_EVAL_NONE); 1160ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "x", ncompx, CEED_EVAL_INTERP); 1161ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ifunction, "v", ncompq, CEED_EVAL_INTERP); 1162ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ifunction, "dv", ncompq*dim, CEED_EVAL_GRAD); 1163ccaff030SJeremy L Thompson } 1164ccaff030SJeremy L Thompson 1165ccaff030SJeremy L Thompson // Create the operator that builds the quadrature data for the NS operator 1166ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup); 1167ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1168ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, 1169ccaff030SJeremy L Thompson basisx, CEED_VECTOR_NONE); 1170ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "qdata", restrictqdi, 1171ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1172ccaff030SJeremy L Thompson 1173ccaff030SJeremy L Thompson // Create the operator that sets the ICs 1174ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 1175ccaff030SJeremy L Thompson CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 1176ccaff030SJeremy L Thompson CeedOperatorSetField(op_ics, "q0", restrictq, 1177ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1178ccaff030SJeremy L Thompson 1179ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 1180ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 1181ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1182ccaff030SJeremy L Thompson 1183ccaff030SJeremy L Thompson if (qf_rhs) { // Create the RHS physics operator 1184ccaff030SJeremy L Thompson CeedOperator op; 1185ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op); 1186ccaff030SJeremy L Thompson CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1187ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1188ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdata", restrictqdi, 1189ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 1190ccaff030SJeremy L Thompson CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1191ccaff030SJeremy L Thompson CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1192ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1193ccaff030SJeremy L Thompson user->op_rhs = op; 1194ccaff030SJeremy L Thompson } 1195ccaff030SJeremy L Thompson 1196ccaff030SJeremy L Thompson if (qf_ifunction) { // Create the IFunction operator 1197ccaff030SJeremy L Thompson CeedOperator op; 1198ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ifunction, NULL, NULL, &op); 1199ccaff030SJeremy L Thompson CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1200ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1201ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 1202ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdata", restrictqdi, 1203ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 1204ccaff030SJeremy L Thompson CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1205ccaff030SJeremy L Thompson CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1206ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1207ccaff030SJeremy L Thompson user->op_ifunction = op; 1208ccaff030SJeremy L Thompson } 1209ccaff030SJeremy L Thompson 1210ccaff030SJeremy L Thompson CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup); 1211ccaff030SJeremy L Thompson CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd}; 1212ccaff030SJeremy L Thompson struct Advection2dContext_ ctxAdvection2d = { 1213ccaff030SJeremy L Thompson .CtauS = CtauS, 1214ccaff030SJeremy L Thompson .strong_form = strong_form, 1215ccaff030SJeremy L Thompson .stabilization = stab, 1216ccaff030SJeremy L Thompson }; 1217ccaff030SJeremy L Thompson switch (problemChoice) { 1218ccaff030SJeremy L Thompson case NS_DENSITY_CURRENT: 1219ccaff030SJeremy L Thompson if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxNS, sizeof ctxNS); 1220ccaff030SJeremy L Thompson if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxNS, 1221ccaff030SJeremy L Thompson sizeof ctxNS); 1222ccaff030SJeremy L Thompson break; 1223ccaff030SJeremy L Thompson case NS_ADVECTION: 1224ccaff030SJeremy L Thompson case NS_ADVECTION2D: 1225ccaff030SJeremy L Thompson if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxAdvection2d, 1226ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 1227ccaff030SJeremy L Thompson if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxAdvection2d, 1228ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 1229ccaff030SJeremy L Thompson } 1230ccaff030SJeremy L Thompson 1231ccaff030SJeremy L Thompson // Set up PETSc context 1232ccaff030SJeremy L Thompson // Set up units structure 1233ccaff030SJeremy L Thompson units->meter = meter; 1234ccaff030SJeremy L Thompson units->kilogram = kilogram; 1235ccaff030SJeremy L Thompson units->second = second; 1236ccaff030SJeremy L Thompson units->Kelvin = Kelvin; 1237ccaff030SJeremy L Thompson units->Pascal = Pascal; 1238ccaff030SJeremy L Thompson units->JperkgK = JperkgK; 1239ccaff030SJeremy L Thompson units->mpersquareds = mpersquareds; 1240ccaff030SJeremy L Thompson units->WpermK = WpermK; 1241ccaff030SJeremy L Thompson units->kgpercubicm = kgpercubicm; 1242ccaff030SJeremy L Thompson units->kgpersquaredms = kgpersquaredms; 1243ccaff030SJeremy L Thompson units->Joulepercubicm = Joulepercubicm; 1244ccaff030SJeremy L Thompson 1245ccaff030SJeremy L Thompson // Set up user structure 1246ccaff030SJeremy L Thompson user->comm = comm; 1247ccaff030SJeremy L Thompson user->outputfreq = outputfreq; 1248ccaff030SJeremy L Thompson user->contsteps = contsteps; 1249ccaff030SJeremy L Thompson user->units = units; 1250ccaff030SJeremy L Thompson user->dm = dm; 1251ccaff030SJeremy L Thompson user->dmviz = dmviz; 1252ccaff030SJeremy L Thompson user->interpviz = interpviz; 1253ccaff030SJeremy L Thompson user->ceed = ceed; 1254ccaff030SJeremy L Thompson 1255ccaff030SJeremy L Thompson // Calculate qdata and ICs 1256ccaff030SJeremy L Thompson // Set up state global and local vectors 1257ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1258ccaff030SJeremy L Thompson 1259ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1260ccaff030SJeremy L Thompson 1261ccaff030SJeremy L Thompson // Apply Setup Ceed Operators 1262ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 1263ccaff030SJeremy L Thompson CeedOperatorApply(op_setup, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 1264ccaff030SJeremy L Thompson ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1265ccaff030SJeremy L Thompson user->M); CHKERRQ(ierr); 1266ccaff030SJeremy L Thompson 126727dd567eSvaleriabarra ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 12683e4faa5aSvaleriabarra &ctxSetup, 0.0); CHKERRQ(ierr); 1269ccaff030SJeremy L Thompson if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1270ccaff030SJeremy L Thompson // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1271ccaff030SJeremy L Thompson // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1272ccaff030SJeremy L Thompson // slower execution. 1273ccaff030SJeremy L Thompson Vec Qbc; 1274ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1275ccaff030SJeremy L Thompson ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1276ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1277ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1278ccaff030SJeremy L Thompson ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1279ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1280ccaff030SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)dm, 1281380e0a58Svaleriabarra "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 1282380e0a58Svaleriabarra CHKERRQ(ierr); 1283ccaff030SJeremy L Thompson } 1284ccaff030SJeremy L Thompson 1285ccaff030SJeremy L Thompson MPI_Comm_rank(comm, &rank); 1286ccaff030SJeremy L Thompson if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);} 1287ccaff030SJeremy L Thompson // Gather initial Q values 1288ccaff030SJeremy L Thompson // In case of continuation of simulation, set up initial values from binary file 1289ccaff030SJeremy L Thompson if (contsteps) { // continue from existent solution 1290ccaff030SJeremy L Thompson PetscViewer viewer; 1291ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1292ccaff030SJeremy L Thompson // Read input 1293ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1294ccaff030SJeremy L Thompson user->outputfolder); 1295ccaff030SJeremy L Thompson CHKERRQ(ierr); 1296ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1297ccaff030SJeremy L Thompson CHKERRQ(ierr); 1298ccaff030SJeremy L Thompson ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1299ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1300ccaff030SJeremy L Thompson } 1301ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1302ccaff030SJeremy L Thompson 1303ccaff030SJeremy L Thompson // Create and setup TS 1304ccaff030SJeremy L Thompson ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1305ccaff030SJeremy L Thompson ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1306ccaff030SJeremy L Thompson if (implicit) { 1307ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1308ccaff030SJeremy L Thompson if (user->op_ifunction) { 1309ccaff030SJeremy L Thompson ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1310ccaff030SJeremy L Thompson } else { // Implicit integrators can fall back to using an RHSFunction 1311ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1312ccaff030SJeremy L Thompson } 1313ccaff030SJeremy L Thompson } else { 1314ccaff030SJeremy L Thompson if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL, 1315ccaff030SJeremy L Thompson "Problem does not provide RHSFunction"); 1316ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1317ccaff030SJeremy L Thompson ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1318ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1319ccaff030SJeremy L Thompson } 1320ccaff030SJeremy L Thompson ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1321ccaff030SJeremy L Thompson ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1322ccaff030SJeremy L Thompson ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 13238ca749c7Svaleriabarra if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);} 1324ccaff030SJeremy L Thompson ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1325ccaff030SJeremy L Thompson ierr = TSAdaptSetStepLimits(adapt, 1326ccaff030SJeremy L Thompson 1.e-12 * units->second, 1327ccaff030SJeremy L Thompson 1.e2 * units->second); CHKERRQ(ierr); 1328ccaff030SJeremy L Thompson ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1329ccaff030SJeremy L Thompson if (!contsteps) { // print initial condition 13308ca749c7Svaleriabarra if (testChoice == TEST_NONE) { 1331ccaff030SJeremy L Thompson ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1332ccaff030SJeremy L Thompson } 1333ccaff030SJeremy L Thompson } else { // continue from time of last output 1334ccaff030SJeremy L Thompson PetscReal time; 1335ccaff030SJeremy L Thompson PetscInt count; 1336ccaff030SJeremy L Thompson PetscViewer viewer; 1337ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1338ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1339ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 1340ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1341ccaff030SJeremy L Thompson CHKERRQ(ierr); 1342ccaff030SJeremy L Thompson ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1343ccaff030SJeremy L Thompson CHKERRQ(ierr); 1344ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1345ccaff030SJeremy L Thompson ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1346ccaff030SJeremy L Thompson } 13478ca749c7Svaleriabarra if (testChoice == TEST_NONE) { 1348ccaff030SJeremy L Thompson ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1349ccaff030SJeremy L Thompson } 1350ccaff030SJeremy L Thompson 1351ccaff030SJeremy L Thompson // Solve 1352ccaff030SJeremy L Thompson start = MPI_Wtime(); 1353ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1354ccaff030SJeremy L Thompson ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1355ccaff030SJeremy L Thompson cpu_time_used = MPI_Wtime() - start; 1356ccaff030SJeremy L Thompson ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1357ccaff030SJeremy L Thompson ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1358ccaff030SJeremy L Thompson comm); CHKERRQ(ierr); 13598ca749c7Svaleriabarra if (testChoice == TEST_NONE) { 1360ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1361380e0a58Svaleriabarra "Time taken for solution (sec): %g\n", 1362ccaff030SJeremy L Thompson (double)cpu_time_used); CHKERRQ(ierr); 1363ccaff030SJeremy L Thompson } 1364ccaff030SJeremy L Thompson 1365ccaff030SJeremy L Thompson // Get error 13668ca749c7Svaleriabarra if (problem->non_zero_time && testChoice == TEST_NONE) { 1367ccaff030SJeremy L Thompson Vec Qexact, Qexactloc; 1368ccaff030SJeremy L Thompson PetscReal norm; 1369ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1370ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1371ccaff030SJeremy L Thompson ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1372ccaff030SJeremy L Thompson 137327dd567eSvaleriabarra ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 1374ccaff030SJeremy L Thompson restrictq, &ctxSetup, ftime); CHKERRQ(ierr); 1375ccaff030SJeremy L Thompson 1376ccaff030SJeremy L Thompson ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1377ccaff030SJeremy L Thompson ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1378ccaff030SJeremy L Thompson CeedVectorDestroy(&q0ceed); 1379ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1380ccaff030SJeremy L Thompson "Max Error: %g\n", 1381ccaff030SJeremy L Thompson (double)norm); CHKERRQ(ierr); 1382357222b3Svaleriabarra // Clean up vectors 1383357222b3Svaleriabarra ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1384357222b3Svaleriabarra ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1385ccaff030SJeremy L Thompson } 1386ccaff030SJeremy L Thompson 1387ccaff030SJeremy L Thompson // Output Statistics 1388ccaff030SJeremy L Thompson ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr); 13898ca749c7Svaleriabarra if (testChoice == TEST_NONE) { 1390ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1391ccaff030SJeremy L Thompson "Time integrator took %D time steps to reach final time %g\n", 1392ccaff030SJeremy L Thompson steps, (double)ftime); CHKERRQ(ierr); 1393ccaff030SJeremy L Thompson } 1394d5424b73Svaleriabarra // Output numerical values from command line 1395d5424b73Svaleriabarra ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr); 1396ccaff030SJeremy L Thompson 13972f9599b8Svaleriabarra // compare reference solution values with current run 13982f9599b8Svaleriabarra if (testChoice != TEST_NONE) { 13999cf88b28Svaleriabarra PetscViewer viewer; 14009cf88b28Svaleriabarra // Read reference file 14019cf88b28Svaleriabarra Vec Qref; 14029cf88b28Svaleriabarra PetscReal error, Qrefnorm; 14033e4faa5aSvaleriabarra ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr); 14048ca749c7Svaleriabarra ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer); 14059cf88b28Svaleriabarra CHKERRQ(ierr); 14069cf88b28Svaleriabarra ierr = VecLoad(Qref, viewer); CHKERRQ(ierr); 14079cf88b28Svaleriabarra ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 14089cf88b28Svaleriabarra 14099cf88b28Svaleriabarra // Compute error with respect to reference solution 14109cf88b28Svaleriabarra ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr); 14119cf88b28Svaleriabarra ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr); 14129cf88b28Svaleriabarra ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr); 14139cf88b28Svaleriabarra ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr); 14149cf88b28Svaleriabarra ierr = VecDestroy(&Qref); CHKERRQ(ierr); 14159cf88b28Svaleriabarra // Check error 14168ca749c7Svaleriabarra if (error > test->testtol) { 14179cf88b28Svaleriabarra ierr = PetscPrintf(PETSC_COMM_WORLD, 14189cf88b28Svaleriabarra "Test failed with error norm %g\n", 14199cf88b28Svaleriabarra (double)error); CHKERRQ(ierr); 14209cf88b28Svaleriabarra } 14213e4faa5aSvaleriabarra ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 14229cf88b28Svaleriabarra } 14239cf88b28Svaleriabarra 1424ccaff030SJeremy L Thompson // Clean up libCEED 1425ccaff030SJeremy L Thompson CeedVectorDestroy(&qdata); 1426ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qceed); 1427ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qdotceed); 1428ccaff030SJeremy L Thompson CeedVectorDestroy(&user->gceed); 1429ccaff030SJeremy L Thompson CeedVectorDestroy(&xcorners); 1430ccaff030SJeremy L Thompson CeedBasisDestroy(&basisq); 1431ccaff030SJeremy L Thompson CeedBasisDestroy(&basisx); 1432ccaff030SJeremy L Thompson CeedBasisDestroy(&basisxc); 1433ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictq); 1434ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictx); 1435ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictqdi); 1436ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_setup); 1437ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ics); 1438ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_rhs); 1439ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ifunction); 1440ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_setup); 1441ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_ics); 1442ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_rhs); 1443ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_ifunction); 1444ccaff030SJeremy L Thompson CeedDestroy(&ceed); 1445ccaff030SJeremy L Thompson 1446ccaff030SJeremy L Thompson // Clean up PETSc 1447ccaff030SJeremy L Thompson ierr = VecDestroy(&Q); CHKERRQ(ierr); 1448ccaff030SJeremy L Thompson ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1449ccaff030SJeremy L Thompson ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1450ccaff030SJeremy L Thompson ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1451ccaff030SJeremy L Thompson ierr = TSDestroy(&ts); CHKERRQ(ierr); 1452ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1453ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 1454ccaff030SJeremy L Thompson ierr = PetscFree(user); CHKERRQ(ierr); 1455ccaff030SJeremy L Thompson return PetscFinalize(); 1456ccaff030SJeremy L Thompson } 1457ccaff030SJeremy L Thompson 1458