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 // 34ff6701fcSJed Brown //TESTARGS -ceed {ceed_resource} -test -degree 1 35ccaff030SJeremy L Thompson 36ccaff030SJeremy L Thompson /// @file 37ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc 38ccaff030SJeremy L Thompson 39ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 40ccaff030SJeremy L Thompson 41ccaff030SJeremy L Thompson #include <petscts.h> 42ccaff030SJeremy L Thompson #include <petscdmplex.h> 43ccaff030SJeremy L Thompson #include <ceed.h> 44ccaff030SJeremy L Thompson #include <stdbool.h> 45ccaff030SJeremy L Thompson #include <petscsys.h> 46ccaff030SJeremy L Thompson #include "common.h" 47ccaff030SJeremy L Thompson #include "advection.h" 48ccaff030SJeremy L Thompson #include "advection2d.h" 49ccaff030SJeremy L Thompson #include "densitycurrent.h" 50ccaff030SJeremy L Thompson 51ccaff030SJeremy L Thompson // Problem Options 52ccaff030SJeremy L Thompson typedef enum { 53ccaff030SJeremy L Thompson NS_DENSITY_CURRENT = 0, 54ccaff030SJeremy L Thompson NS_ADVECTION = 1, 55ccaff030SJeremy L Thompson NS_ADVECTION2D = 2, 56ccaff030SJeremy L Thompson } problemType; 57ccaff030SJeremy L Thompson static const char *const problemTypes[] = { 58ccaff030SJeremy L Thompson "density_current", 59ccaff030SJeremy L Thompson "advection", 60ccaff030SJeremy L Thompson "advection2d", 61ccaff030SJeremy L Thompson "problemType","NS_",0 62ccaff030SJeremy L Thompson }; 63ccaff030SJeremy L Thompson 64ccaff030SJeremy L Thompson typedef enum { 65ccaff030SJeremy L Thompson STAB_NONE = 0, 66ccaff030SJeremy L Thompson STAB_SU = 1, // Streamline Upwind 67ccaff030SJeremy L Thompson STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 68ccaff030SJeremy L Thompson } StabilizationType; 69ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = { 70ccaff030SJeremy L Thompson "NONE", 71ccaff030SJeremy L Thompson "SU", 72ccaff030SJeremy L Thompson "SUPG", 73ccaff030SJeremy L Thompson "StabilizationType", "STAB_", NULL 74ccaff030SJeremy L Thompson }; 75ccaff030SJeremy L Thompson 76ccaff030SJeremy L Thompson // Problem specific data 77ccaff030SJeremy L Thompson typedef struct { 78ccaff030SJeremy L Thompson CeedInt dim, qdatasize; 79ccaff030SJeremy L Thompson CeedQFunctionUser setup, ics, apply_rhs, apply_ifunction; 80ccaff030SJeremy L Thompson PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 81ccaff030SJeremy L Thompson PetscScalar[], void *); 82ccaff030SJeremy L Thompson const char *setup_loc, *ics_loc, *apply_rhs_loc, *apply_ifunction_loc; 83ccaff030SJeremy L Thompson const bool non_zero_time; 84ccaff030SJeremy L Thompson } problemData; 85ccaff030SJeremy L Thompson 86ccaff030SJeremy L Thompson problemData problemOptions[] = { 87ccaff030SJeremy L Thompson [NS_DENSITY_CURRENT] = { 88ccaff030SJeremy L Thompson .dim = 3, 89ccaff030SJeremy L Thompson .qdatasize = 10, 90ccaff030SJeremy L Thompson .setup = Setup, 91ccaff030SJeremy L Thompson .setup_loc = Setup_loc, 92ccaff030SJeremy L Thompson .ics = ICsDC, 93ccaff030SJeremy L Thompson .ics_loc = ICsDC_loc, 94ccaff030SJeremy L Thompson .apply_rhs = DC, 95ccaff030SJeremy L Thompson .apply_rhs_loc = DC_loc, 96ccaff030SJeremy L Thompson .apply_ifunction = IFunction_DC, 97ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_DC_loc, 98ccaff030SJeremy L Thompson .bc = Exact_DC, 99ccaff030SJeremy L Thompson .non_zero_time = false, 100ccaff030SJeremy L Thompson }, 101ccaff030SJeremy L Thompson [NS_ADVECTION] = { 102ccaff030SJeremy L Thompson .dim = 3, 103ccaff030SJeremy L Thompson .qdatasize = 10, 104ccaff030SJeremy L Thompson .setup = Setup, 105ccaff030SJeremy L Thompson .setup_loc = Setup_loc, 106ccaff030SJeremy L Thompson .ics = ICsAdvection, 107ccaff030SJeremy L Thompson .ics_loc = ICsAdvection_loc, 108ccaff030SJeremy L Thompson .apply_rhs = Advection, 109ccaff030SJeremy L Thompson .apply_rhs_loc = Advection_loc, 110ccaff030SJeremy L Thompson .apply_ifunction = IFunction_Advection, 111ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_Advection_loc, 112ccaff030SJeremy L Thompson .bc = Exact_Advection, 113ec9acdb8Svaleriabarra .non_zero_time = false, 114ccaff030SJeremy L Thompson }, 115ccaff030SJeremy L Thompson [NS_ADVECTION2D] = { 116ccaff030SJeremy L Thompson .dim = 2, 117ccaff030SJeremy L Thompson .qdatasize = 5, 118ccaff030SJeremy L Thompson .setup = Setup2d, 119ccaff030SJeremy L Thompson .setup_loc = Setup2d_loc, 120ccaff030SJeremy L Thompson .ics = ICsAdvection2d, 121ccaff030SJeremy L Thompson .ics_loc = ICsAdvection2d_loc, 122ccaff030SJeremy L Thompson .apply_rhs = Advection2d, 123ccaff030SJeremy L Thompson .apply_rhs_loc = Advection2d_loc, 124ccaff030SJeremy L Thompson .apply_ifunction = IFunction_Advection2d, 125ccaff030SJeremy L Thompson .apply_ifunction_loc = IFunction_Advection2d_loc, 126ccaff030SJeremy L Thompson .bc = Exact_Advection2d, 127ccaff030SJeremy L Thompson .non_zero_time = true, 128ccaff030SJeremy L Thompson }, 129ccaff030SJeremy L Thompson }; 130ccaff030SJeremy L Thompson 131ccaff030SJeremy L Thompson // PETSc user data 132ccaff030SJeremy L Thompson typedef struct User_ *User; 133ccaff030SJeremy L Thompson typedef struct Units_ *Units; 134ccaff030SJeremy L Thompson 135ccaff030SJeremy L Thompson struct User_ { 136ccaff030SJeremy L Thompson MPI_Comm comm; 137ccaff030SJeremy L Thompson PetscInt outputfreq; 138ccaff030SJeremy L Thompson DM dm; 139ccaff030SJeremy L Thompson DM dmviz; 140ccaff030SJeremy L Thompson Mat interpviz; 141ccaff030SJeremy L Thompson Ceed ceed; 142ccaff030SJeremy L Thompson Units units; 143ccaff030SJeremy L Thompson CeedVector qceed, qdotceed, gceed; 144ccaff030SJeremy L Thompson CeedOperator op_rhs, op_ifunction; 145ccaff030SJeremy L Thompson Vec M; 146ccaff030SJeremy L Thompson char outputfolder[PETSC_MAX_PATH_LEN]; 147ccaff030SJeremy L Thompson PetscInt contsteps; 148ccaff030SJeremy L Thompson }; 149ccaff030SJeremy L Thompson 150ccaff030SJeremy L Thompson struct Units_ { 151ccaff030SJeremy L Thompson // fundamental units 152ccaff030SJeremy L Thompson PetscScalar meter; 153ccaff030SJeremy L Thompson PetscScalar kilogram; 154ccaff030SJeremy L Thompson PetscScalar second; 155ccaff030SJeremy L Thompson PetscScalar Kelvin; 156ccaff030SJeremy L Thompson // derived units 157ccaff030SJeremy L Thompson PetscScalar Pascal; 158ccaff030SJeremy L Thompson PetscScalar JperkgK; 159ccaff030SJeremy L Thompson PetscScalar mpersquareds; 160ccaff030SJeremy L Thompson PetscScalar WpermK; 161ccaff030SJeremy L Thompson PetscScalar kgpercubicm; 162ccaff030SJeremy L Thompson PetscScalar kgpersquaredms; 163ccaff030SJeremy L Thompson PetscScalar Joulepercubicm; 164ccaff030SJeremy L Thompson }; 165ccaff030SJeremy L Thompson 166ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC; 167ccaff030SJeremy L Thompson struct SimpleBC_ { 168ccaff030SJeremy L Thompson PetscInt nwall, nslip[3]; 169ccaff030SJeremy L Thompson PetscInt walls[10], slips[3][10]; 170ccaff030SJeremy L Thompson }; 171ccaff030SJeremy L Thompson 172ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1). 173ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) { 174ccaff030SJeremy L Thompson return i >= 0 ? i : -(i+1); 175ccaff030SJeremy L Thompson } 176ccaff030SJeremy L Thompson 177ccaff030SJeremy L Thompson // Utility function to create local CEED restriction 178ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 179ccaff030SJeremy L Thompson CeedElemRestriction *Erestrict) { 180ccaff030SJeremy L Thompson 181ccaff030SJeremy L Thompson PetscSection section; 182ccaff030SJeremy L Thompson PetscInt c, cStart, cEnd, Nelem, Ndof, *erestrict, eoffset, nfields, dim; 183ccaff030SJeremy L Thompson PetscErrorCode ierr; 184ccaff030SJeremy L Thompson Vec Uloc; 185ccaff030SJeremy L Thompson 186ccaff030SJeremy L Thompson PetscFunctionBeginUser; 187ccaff030SJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 188ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm,§ion); CHKERRQ(ierr); 189ccaff030SJeremy L Thompson ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 190ccaff030SJeremy L Thompson PetscInt ncomp[nfields], fieldoff[nfields+1]; 191ccaff030SJeremy L Thompson fieldoff[0] = 0; 192ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 193ccaff030SJeremy L Thompson ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 194ccaff030SJeremy L Thompson fieldoff[f+1] = fieldoff[f] + ncomp[f]; 195ccaff030SJeremy L Thompson } 196ccaff030SJeremy L Thompson 197ccaff030SJeremy L Thompson ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 198ccaff030SJeremy L Thompson Nelem = cEnd - cStart; 199ccaff030SJeremy L Thompson ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 200ccaff030SJeremy L Thompson for (c=cStart,eoffset=0; c<cEnd; c++) { 201ccaff030SJeremy L Thompson PetscInt numindices, *indices, nnodes; 202ccaff030SJeremy L Thompson ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices, 203ccaff030SJeremy L Thompson &indices, NULL); CHKERRQ(ierr); 204ccaff030SJeremy L Thompson if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 205ccaff030SJeremy L Thompson PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 206ccaff030SJeremy L Thompson c); 207ccaff030SJeremy L Thompson nnodes = numindices / fieldoff[nfields]; 208ccaff030SJeremy L Thompson for (PetscInt i=0; i<nnodes; i++) { 209ccaff030SJeremy L Thompson // Check that indices are blocked by node and thus can be coalesced as a single field with 210ccaff030SJeremy L Thompson // fieldoff[nfields] = sum(ncomp) components. 211ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 212ccaff030SJeremy L Thompson for (PetscInt j=0; j<ncomp[f]; j++) { 213ccaff030SJeremy L Thompson if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j]) 214ccaff030SJeremy L Thompson != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j) 215ccaff030SJeremy L Thompson SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 216ccaff030SJeremy L Thompson "Cell %D closure indices not interlaced for node %D field %D component %D", 217ccaff030SJeremy L Thompson c, i, f, j); 218ccaff030SJeremy L Thompson } 219ccaff030SJeremy L Thompson } 220ccaff030SJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 221ccaff030SJeremy L Thompson PetscInt loc = Involute(indices[i*ncomp[0]]); 222ccaff030SJeremy L Thompson erestrict[eoffset++] = loc / fieldoff[nfields]; 223ccaff030SJeremy L Thompson } 224ccaff030SJeremy L Thompson ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices, 225ccaff030SJeremy L Thompson &indices, NULL); CHKERRQ(ierr); 226ccaff030SJeremy L Thompson } 227ccaff030SJeremy L Thompson if (eoffset != Nelem*PetscPowInt(P, dim)) SETERRQ3(PETSC_COMM_SELF, 228ccaff030SJeremy L Thompson PETSC_ERR_LIB, "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 229ccaff030SJeremy L Thompson PetscPowInt(P, dim),eoffset); 230ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 231ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 232ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 233ccaff030SJeremy L Thompson CeedElemRestrictionCreate(ceed, CEED_INTERLACED, Nelem, PetscPowInt(P, dim), 234ccaff030SJeremy L Thompson Ndof/fieldoff[nfields], fieldoff[nfields], 235ccaff030SJeremy L Thompson CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, Erestrict); 236ccaff030SJeremy L Thompson ierr = PetscFree(erestrict); CHKERRQ(ierr); 237ccaff030SJeremy L Thompson PetscFunctionReturn(0); 238ccaff030SJeremy L Thompson } 239ccaff030SJeremy L Thompson 240ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 241ccaff030SJeremy L Thompson PetscErrorCode ierr; 242ccaff030SJeremy L Thompson PetscInt m; 243ccaff030SJeremy L Thompson 244ccaff030SJeremy L Thompson PetscFunctionBeginUser; 245ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 246ccaff030SJeremy L Thompson ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 247ccaff030SJeremy L Thompson PetscFunctionReturn(0); 248ccaff030SJeremy L Thompson } 249ccaff030SJeremy L Thompson 250ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) { 251ccaff030SJeremy L Thompson PetscErrorCode ierr; 252ccaff030SJeremy L Thompson PetscInt mceed,mpetsc; 253ccaff030SJeremy L Thompson PetscScalar *a; 254ccaff030SJeremy L Thompson 255ccaff030SJeremy L Thompson PetscFunctionBeginUser; 256ccaff030SJeremy L Thompson ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 257ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 258ccaff030SJeremy L Thompson if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 259380e0a58Svaleriabarra "Cannot place PETSc Vec of length %D in CeedVector of length %D", 260380e0a58Svaleriabarra mpetsc, mceed); 261ccaff030SJeremy L Thompson ierr = VecGetArray(p, &a); CHKERRQ(ierr); 262ccaff030SJeremy L Thompson CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 263ccaff030SJeremy L Thompson PetscFunctionReturn(0); 264ccaff030SJeremy L Thompson } 265ccaff030SJeremy L Thompson 266ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 267ccaff030SJeremy L Thompson PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 268ccaff030SJeremy L Thompson Vec cellGeomFVM, Vec gradFVM) { 269ccaff030SJeremy L Thompson PetscErrorCode ierr; 270ccaff030SJeremy L Thompson Vec Qbc; 271ccaff030SJeremy L Thompson 272ccaff030SJeremy L Thompson PetscFunctionBegin; 273ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 274ccaff030SJeremy L Thompson ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 275ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 276ccaff030SJeremy L Thompson PetscFunctionReturn(0); 277ccaff030SJeremy L Thompson } 278ccaff030SJeremy L Thompson 279ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u) 280ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G 281ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 282ccaff030SJeremy L Thompson PetscErrorCode ierr; 283ccaff030SJeremy L Thompson User user = *(User *)userData; 284ccaff030SJeremy L Thompson PetscScalar *q, *g; 285ccaff030SJeremy L Thompson Vec Qloc, Gloc; 286ccaff030SJeremy L Thompson 287ccaff030SJeremy L Thompson // Global-to-local 288ccaff030SJeremy L Thompson PetscFunctionBeginUser; 289ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 290ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 291ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 292ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 293ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 294ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 295ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 296ccaff030SJeremy L Thompson 297ccaff030SJeremy L Thompson // Ceed Vectors 298ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 299ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 300ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 301ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 302ccaff030SJeremy L Thompson 303ccaff030SJeremy L Thompson // Apply CEED operator 304ccaff030SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 305ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 306ccaff030SJeremy L Thompson 307ccaff030SJeremy L Thompson // Restore vectors 308ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 309ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 310ccaff030SJeremy L Thompson 311ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 312ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 313ccaff030SJeremy L Thompson 314ccaff030SJeremy L Thompson // Inverse of the lumped mass matrix 315ccaff030SJeremy L Thompson ierr = VecPointwiseMult(G, G, user->M); // M is Minv 316ccaff030SJeremy L Thompson CHKERRQ(ierr); 317ccaff030SJeremy L Thompson 318ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 319ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 320ccaff030SJeremy L Thompson PetscFunctionReturn(0); 321ccaff030SJeremy L Thompson } 322ccaff030SJeremy L Thompson 323ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 324ccaff030SJeremy L Thompson void *userData) { 325ccaff030SJeremy L Thompson PetscErrorCode ierr; 326ccaff030SJeremy L Thompson User user = *(User *)userData; 327ccaff030SJeremy L Thompson const PetscScalar *q, *qdot; 328ccaff030SJeremy L Thompson PetscScalar *g; 329ccaff030SJeremy L Thompson Vec Qloc, Qdotloc, Gloc; 330ccaff030SJeremy L Thompson 331ccaff030SJeremy L Thompson // Global-to-local 332ccaff030SJeremy L Thompson PetscFunctionBeginUser; 333ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 334ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 335ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 336ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 337ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 338ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 339ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 340ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 341ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 342ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 343ccaff030SJeremy L Thompson 344ccaff030SJeremy L Thompson // Ceed Vectors 345ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 346ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 347ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 348ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 349ccaff030SJeremy L Thompson (PetscScalar *)q); 350ccaff030SJeremy L Thompson CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 351ccaff030SJeremy L Thompson (PetscScalar *)qdot); 352ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 353ccaff030SJeremy L Thompson 354ccaff030SJeremy L Thompson // Apply CEED operator 355ccaff030SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 356ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 357ccaff030SJeremy L Thompson 358ccaff030SJeremy L Thompson // Restore vectors 359ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 360ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 361ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 362ccaff030SJeremy L Thompson 363ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 364ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 365ccaff030SJeremy L Thompson 366ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 367ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 368ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 369ccaff030SJeremy L Thompson PetscFunctionReturn(0); 370ccaff030SJeremy L Thompson } 371ccaff030SJeremy L Thompson 372ccaff030SJeremy L Thompson // User provided TS Monitor 373ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 374ccaff030SJeremy L Thompson Vec Q, void *ctx) { 375ccaff030SJeremy L Thompson User user = ctx; 376ccaff030SJeremy L Thompson Vec Qloc; 377ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 378ccaff030SJeremy L Thompson PetscViewer viewer; 379ccaff030SJeremy L Thompson PetscErrorCode ierr; 380ccaff030SJeremy L Thompson 381ccaff030SJeremy L Thompson // Set up output 382ccaff030SJeremy L Thompson PetscFunctionBeginUser; 383ccaff030SJeremy L Thompson // Print every 'outputfreq' steps 384ccaff030SJeremy L Thompson if (stepno % user->outputfreq != 0) 385ccaff030SJeremy L Thompson PetscFunctionReturn(0); 386ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 387ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 388ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 389ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 390ccaff030SJeremy L Thompson 391ccaff030SJeremy L Thompson // Output 392ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 393ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 394ccaff030SJeremy L Thompson CHKERRQ(ierr); 395ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 396ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 397ccaff030SJeremy L Thompson ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 398ccaff030SJeremy L Thompson if (user->dmviz) { 399ccaff030SJeremy L Thompson Vec Qrefined, Qrefined_loc; 400ccaff030SJeremy L Thompson char filepath_refined[PETSC_MAX_PATH_LEN]; 401ccaff030SJeremy L Thompson PetscViewer viewer_refined; 402ccaff030SJeremy L Thompson 403ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 404ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 405ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 406ccaff030SJeremy L Thompson CHKERRQ(ierr); 407ccaff030SJeremy L Thompson ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 408ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 409ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 410ccaff030SJeremy L Thompson CHKERRQ(ierr); 411ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 412ccaff030SJeremy L Thompson "%s/nsrefined-%03D.vtu", 413ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 414ccaff030SJeremy L Thompson CHKERRQ(ierr); 415ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 416ccaff030SJeremy L Thompson filepath_refined, 417ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 418ccaff030SJeremy L Thompson ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 419ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 420ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 421ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 422ccaff030SJeremy L Thompson } 423ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 424ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 425ccaff030SJeremy L Thompson 426ccaff030SJeremy L Thompson // Save data in a binary file for continuation of simulations 427ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 428ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 429ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 430ccaff030SJeremy L Thompson CHKERRQ(ierr); 431ccaff030SJeremy L Thompson ierr = VecView(Q, viewer); CHKERRQ(ierr); 432ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 433ccaff030SJeremy L Thompson 434ccaff030SJeremy L Thompson // Save time stamp 435ccaff030SJeremy L Thompson // Dimensionalize time back 436ccaff030SJeremy L Thompson time /= user->units->second; 437ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 438ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 439ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 440ccaff030SJeremy L Thompson CHKERRQ(ierr); 441ccaff030SJeremy L Thompson #if PETSC_VERSION_GE(3,13,0) 442ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 443ccaff030SJeremy L Thompson #else 444ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 445ccaff030SJeremy L Thompson #endif 446ccaff030SJeremy L Thompson CHKERRQ(ierr); 447ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 448ccaff030SJeremy L Thompson 449ccaff030SJeremy L Thompson PetscFunctionReturn(0); 450ccaff030SJeremy L Thompson } 451ccaff030SJeremy L Thompson 45227dd567eSvaleriabarra static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 453ccaff030SJeremy L Thompson CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 454ccaff030SJeremy L Thompson CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) { 455ccaff030SJeremy L Thompson PetscErrorCode ierr; 456ccaff030SJeremy L Thompson CeedVector multlvec; 457ccaff030SJeremy L Thompson Vec Multiplicity, MultiplicityLoc; 458ccaff030SJeremy L Thompson 459ccaff030SJeremy L Thompson ctxSetup->time = time; 460ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 461ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 462ccaff030SJeremy L Thompson CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 463ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 464ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 465ccaff030SJeremy L Thompson 466ccaff030SJeremy L Thompson // Fix multiplicity for output of ICs 467ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 468ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 469ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 470ccaff030SJeremy L Thompson CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 471ccaff030SJeremy L Thompson CeedVectorDestroy(&multlvec); 472ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 473ccaff030SJeremy L Thompson ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 474ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 475ccaff030SJeremy L Thompson CHKERRQ(ierr); 476ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 477ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 478ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 479ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 480ccaff030SJeremy L Thompson 481ccaff030SJeremy L Thompson PetscFunctionReturn(0); 482ccaff030SJeremy L Thompson } 483ccaff030SJeremy L Thompson 484ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 485ccaff030SJeremy L Thompson CeedElemRestriction restrictq, CeedBasis basisq, 486ccaff030SJeremy L Thompson CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 487ccaff030SJeremy L Thompson PetscErrorCode ierr; 488ccaff030SJeremy L Thompson CeedQFunction qf_mass; 489ccaff030SJeremy L Thompson CeedOperator op_mass; 490ccaff030SJeremy L Thompson CeedVector mceed; 491ccaff030SJeremy L Thompson Vec Mloc; 492ccaff030SJeremy L Thompson CeedInt ncompq, qdatasize; 493ccaff030SJeremy L Thompson 494ccaff030SJeremy L Thompson PetscFunctionBeginUser; 495ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 496ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 497ccaff030SJeremy L Thompson // Create the Q-function that defines the action of the mass operator 498ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 499ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 500ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 501ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 502ccaff030SJeremy L Thompson 503ccaff030SJeremy L Thompson // Create the mass operator 504ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 505ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 506ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", restrictqdi, 507ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 508ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 509ccaff030SJeremy L Thompson 510ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 511ccaff030SJeremy L Thompson ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 512ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 513ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 514ccaff030SJeremy L Thompson 515ccaff030SJeremy L Thompson { 516ccaff030SJeremy L Thompson // Compute a lumped mass matrix 517ccaff030SJeremy L Thompson CeedVector onesvec; 518ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 519ccaff030SJeremy L Thompson CeedVectorSetValue(onesvec, 1.0); 520ccaff030SJeremy L Thompson CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 521ccaff030SJeremy L Thompson CeedVectorDestroy(&onesvec); 522ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_mass); 523ccaff030SJeremy L Thompson CeedVectorDestroy(&mceed); 524ccaff030SJeremy L Thompson } 525ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_mass); 526ccaff030SJeremy L Thompson 527ccaff030SJeremy L Thompson ierr = VecZeroEntries(M); CHKERRQ(ierr); 528ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 529ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 530ccaff030SJeremy L Thompson 531ccaff030SJeremy L Thompson // Invert diagonally lumped mass vector for RHS function 532ccaff030SJeremy L Thompson ierr = VecReciprocal(M); CHKERRQ(ierr); 533ccaff030SJeremy L Thompson PetscFunctionReturn(0); 534ccaff030SJeremy L Thompson } 535ccaff030SJeremy L Thompson 536ff6701fcSJed Brown PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 537ff6701fcSJed Brown SimpleBC bc, void *ctxSetup) { 538ccaff030SJeremy L Thompson PetscErrorCode ierr; 539ccaff030SJeremy L Thompson 540ccaff030SJeremy L Thompson PetscFunctionBeginUser; 541ccaff030SJeremy L Thompson { 542ccaff030SJeremy L Thompson // Configure the finite element space and boundary conditions 543ccaff030SJeremy L Thompson PetscFE fe; 544ccaff030SJeremy L Thompson PetscSpace fespace; 545ccaff030SJeremy L Thompson PetscInt ncompq = 5; 546ff6701fcSJed Brown ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 547ff6701fcSJed Brown PETSC_FALSE, degree, PETSC_DECIDE, 548ff6701fcSJed Brown &fe); 549ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 550ccaff030SJeremy L Thompson ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr); 551ccaff030SJeremy L Thompson ierr = DMCreateDS(dm); CHKERRQ(ierr); 552*35d2ed1bSvaleriabarra // Wall boundary conditions are zero energy density and zero flux for 553*35d2ed1bSvaleriabarra // velocity in advection/advection2d, and zero velocity and zero flux 554*35d2ed1bSvaleriabarra // for mass density and energy density in density_current 555ccaff030SJeremy L Thompson { 556*35d2ed1bSvaleriabarra if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 557*35d2ed1bSvaleriabarra PetscInt comps[1] = {4}; 558*35d2ed1bSvaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 559*35d2ed1bSvaleriabarra 1, comps, (void(*)(void))problem->bc, 560*35d2ed1bSvaleriabarra bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 561*35d2ed1bSvaleriabarra } else if (problem->bc == Exact_DC) { 562ccaff030SJeremy L Thompson PetscInt comps[3] = {1, 2, 3}; 563ccaff030SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 564ccaff030SJeremy L Thompson 3, comps, (void(*)(void))problem->bc, 565ccaff030SJeremy L Thompson bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 566ccaff030SJeremy L Thompson } 567*35d2ed1bSvaleriabarra } 568ccaff030SJeremy L Thompson { 569ccaff030SJeremy L Thompson PetscInt comps[1] = {1}; 570ccaff030SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 571ccaff030SJeremy L Thompson 1, comps, (void(*)(void))NULL, bc->nslip[0], 572ccaff030SJeremy L Thompson bc->slips[0], ctxSetup); CHKERRQ(ierr); 573ccaff030SJeremy L Thompson comps[0] = 2; 574ccaff030SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 575ccaff030SJeremy L Thompson 1, comps, (void(*)(void))NULL, bc->nslip[1], 576ccaff030SJeremy L Thompson bc->slips[1], ctxSetup); CHKERRQ(ierr); 577ccaff030SJeremy L Thompson comps[0] = 3; 578ccaff030SJeremy L Thompson ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 579ccaff030SJeremy L Thompson 1, comps, (void(*)(void))NULL, bc->nslip[2], 580ccaff030SJeremy L Thompson bc->slips[2], ctxSetup); CHKERRQ(ierr); 581ccaff030SJeremy L Thompson } 582ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm,PETSC_DETERMINE,NULL); 583ccaff030SJeremy L Thompson CHKERRQ(ierr); 584ccaff030SJeremy L Thompson ierr = PetscFEGetBasisSpace(fe, &fespace); CHKERRQ(ierr); 585ccaff030SJeremy L Thompson ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 586ccaff030SJeremy L Thompson } 587ccaff030SJeremy L Thompson { 588ccaff030SJeremy L Thompson // Empty name for conserved field (because there is only one field) 589ccaff030SJeremy L Thompson PetscSection section; 590ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 591ccaff030SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 592ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 593ccaff030SJeremy L Thompson CHKERRQ(ierr); 594ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 595ccaff030SJeremy L Thompson CHKERRQ(ierr); 596ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 597ccaff030SJeremy L Thompson CHKERRQ(ierr); 598ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 599ccaff030SJeremy L Thompson CHKERRQ(ierr); 600ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 601ccaff030SJeremy L Thompson CHKERRQ(ierr); 602ccaff030SJeremy L Thompson } 603ccaff030SJeremy L Thompson PetscFunctionReturn(0); 604ccaff030SJeremy L Thompson } 605ccaff030SJeremy L Thompson 606ccaff030SJeremy L Thompson int main(int argc, char **argv) { 607ccaff030SJeremy L Thompson PetscInt ierr; 608ccaff030SJeremy L Thompson MPI_Comm comm; 609380e0a58Svaleriabarra DM dm, dmcoord, dmviz; 610ccaff030SJeremy L Thompson Mat interpviz; 611ccaff030SJeremy L Thompson TS ts; 612ccaff030SJeremy L Thompson TSAdapt adapt; 613ccaff030SJeremy L Thompson User user; 614ccaff030SJeremy L Thompson Units units; 615ccaff030SJeremy L Thompson char ceedresource[4096] = "/cpu/self"; 616ccaff030SJeremy L Thompson PetscInt cStart, cEnd, localNelem, lnodes, steps; 617ccaff030SJeremy L Thompson const PetscInt ncompq = 5; 618ccaff030SJeremy L Thompson PetscMPIInt rank; 619ccaff030SJeremy L Thompson PetscScalar ftime; 620ccaff030SJeremy L Thompson Vec Q, Qloc, Xloc; 621ccaff030SJeremy L Thompson Ceed ceed; 622ccaff030SJeremy L Thompson CeedInt numP, numQ; 623ccaff030SJeremy L Thompson CeedVector xcorners, qdata, q0ceed; 624ccaff030SJeremy L Thompson CeedBasis basisx, basisxc, basisq; 625ccaff030SJeremy L Thompson CeedElemRestriction restrictx, restrictxcoord, restrictq, restrictqdi; 626ccaff030SJeremy L Thompson CeedQFunction qf_setup, qf_ics, qf_rhs, qf_ifunction; 627ccaff030SJeremy L Thompson CeedOperator op_setup, op_ics; 628ccaff030SJeremy L Thompson CeedScalar Rd; 629ccaff030SJeremy L Thompson PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 630ccaff030SJeremy L Thompson kgpersquaredms, Joulepercubicm; 631ccaff030SJeremy L Thompson problemType problemChoice; 632ccaff030SJeremy L Thompson problemData *problem = NULL; 633ccaff030SJeremy L Thompson StabilizationType stab; 634cb3e2689Svaleriabarra PetscBool test, implicit; 635cb3e2689Svaleriabarra PetscInt viz_refine = 0; 636ccaff030SJeremy L Thompson struct SimpleBC_ bc = { 637ccaff030SJeremy L Thompson .nwall = 6, 638ccaff030SJeremy L Thompson .walls = {1,2,3,4,5,6}, 639ccaff030SJeremy L Thompson }; 640ccaff030SJeremy L Thompson double start, cpu_time_used; 641ccaff030SJeremy L Thompson 642ccaff030SJeremy L Thompson // Create the libCEED contexts 643ccaff030SJeremy L Thompson PetscScalar meter = 1e-2; // 1 meter in scaled length units 644ccaff030SJeremy L Thompson PetscScalar second = 1e-2; // 1 second in scaled time units 645ccaff030SJeremy L Thompson PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 646ccaff030SJeremy L Thompson PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 647ccaff030SJeremy L Thompson CeedScalar theta0 = 300.; // K 648ccaff030SJeremy L Thompson CeedScalar thetaC = -15.; // K 649ccaff030SJeremy L Thompson CeedScalar P0 = 1.e5; // Pa 650ccaff030SJeremy L Thompson CeedScalar N = 0.01; // 1/s 651ccaff030SJeremy L Thompson CeedScalar cv = 717.; // J/(kg K) 652ccaff030SJeremy L Thompson CeedScalar cp = 1004.; // J/(kg K) 653ccaff030SJeremy L Thompson CeedScalar g = 9.81; // m/s^2 654ccaff030SJeremy L Thompson CeedScalar lambda = -2./3.; // - 655ccaff030SJeremy L Thompson CeedScalar mu = 75.; // Pa s, dynamic viscosity 656ccaff030SJeremy L Thompson // mu = 75 is not physical for air, but is good for numerical stability 657ccaff030SJeremy L Thompson CeedScalar k = 0.02638; // W/(m K) 658ccaff030SJeremy L Thompson CeedScalar CtauS = 0.; // dimensionless 659ccaff030SJeremy L Thompson CeedScalar strong_form = 0.; // [0,1] 660ccaff030SJeremy L Thompson PetscScalar lx = 8000.; // m 661ccaff030SJeremy L Thompson PetscScalar ly = 8000.; // m 662ccaff030SJeremy L Thompson PetscScalar lz = 4000.; // m 663ccaff030SJeremy L Thompson CeedScalar rc = 1000.; // m (Radius of bubble) 664ccaff030SJeremy L Thompson PetscScalar resx = 1000.; // m (resolution in x) 665ccaff030SJeremy L Thompson PetscScalar resy = 1000.; // m (resolution in y) 666ccaff030SJeremy L Thompson PetscScalar resz = 1000.; // m (resolution in z) 667ccaff030SJeremy L Thompson PetscInt outputfreq = 10; // - 668ccaff030SJeremy L Thompson PetscInt contsteps = 0; // - 669ff6701fcSJed Brown PetscInt degree = 1; // - 670ccaff030SJeremy L Thompson PetscInt qextra = 2; // - 671ccaff030SJeremy L Thompson DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, 672ccaff030SJeremy L Thompson DM_BOUNDARY_NONE 673ccaff030SJeremy L Thompson }; 674ccaff030SJeremy L Thompson PetscReal center[3], dc_axis[3] = {0, 0, 0}; 675ccaff030SJeremy L Thompson 676ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 677ccaff030SJeremy L Thompson if (ierr) return ierr; 678ccaff030SJeremy L Thompson 679ccaff030SJeremy L Thompson // Allocate PETSc context 680ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 681ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 682ccaff030SJeremy L Thompson 683ccaff030SJeremy L Thompson // Parse command line options 684ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 685ccaff030SJeremy L Thompson ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 686ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 687ccaff030SJeremy L Thompson ierr = PetscOptionsString("-ceed", "CEED resource specifier", 688ccaff030SJeremy L Thompson NULL, ceedresource, ceedresource, 689ccaff030SJeremy L Thompson sizeof(ceedresource), NULL); CHKERRQ(ierr); 690ccaff030SJeremy L Thompson ierr = PetscOptionsBool("-test", "Run in test mode", 691ccaff030SJeremy L Thompson NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr); 692ccaff030SJeremy L Thompson problemChoice = NS_DENSITY_CURRENT; 693ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 694ccaff030SJeremy L Thompson problemTypes, (PetscEnum)problemChoice, 695ccaff030SJeremy L Thompson (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 696ccaff030SJeremy L Thompson problem = &problemOptions[problemChoice]; 697ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 698ccaff030SJeremy L Thompson StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 699ccaff030SJeremy L Thompson (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 700ccaff030SJeremy L Thompson ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 701ccaff030SJeremy L Thompson NULL, implicit=PETSC_FALSE, &implicit, NULL); 702ccaff030SJeremy L Thompson CHKERRQ(ierr); 7038c662231Svaleriabarra if (!implicit && stab != STAB_NONE) { 7048c662231Svaleriabarra ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 7058c662231Svaleriabarra CHKERRQ(ierr); 7068c662231Svaleriabarra } 707ccaff030SJeremy L Thompson { 708ccaff030SJeremy L Thompson PetscInt len; 709ccaff030SJeremy L Thompson PetscBool flg; 710ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray("-bc_wall", 711ccaff030SJeremy L Thompson "Use wall boundary conditions on this list of faces", 712ccaff030SJeremy L Thompson NULL, bc.walls, 713ccaff030SJeremy L Thompson (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 714ccaff030SJeremy L Thompson &len), &flg); CHKERRQ(ierr); 715ccaff030SJeremy L Thompson if (flg) bc.nwall = len; 716ccaff030SJeremy L Thompson for (PetscInt j=0; j<3; j++) { 717ccaff030SJeremy L Thompson const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 718ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray(flags[j], 719ccaff030SJeremy L Thompson "Use slip boundary conditions on this list of faces", 720ccaff030SJeremy L Thompson NULL, bc.slips[j], 721ccaff030SJeremy L Thompson (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 722ccaff030SJeremy L Thompson &len), &flg); 723ccaff030SJeremy L Thompson CHKERRQ(ierr); 724ccaff030SJeremy L Thompson if (flg) bc.nslip[j] = len; 725ccaff030SJeremy L Thompson } 726ccaff030SJeremy L Thompson } 727cb3e2689Svaleriabarra ierr = PetscOptionsInt("-viz_refine", 728cb3e2689Svaleriabarra "Regular refinement levels for visualization", 729cb3e2689Svaleriabarra NULL, viz_refine, &viz_refine, NULL); 730ccaff030SJeremy L Thompson CHKERRQ(ierr); 731ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 732ccaff030SJeremy L Thompson NULL, meter, &meter, NULL); CHKERRQ(ierr); 733ccaff030SJeremy L Thompson meter = fabs(meter); 734ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 735ccaff030SJeremy L Thompson NULL, second, &second, NULL); CHKERRQ(ierr); 736ccaff030SJeremy L Thompson second = fabs(second); 737ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 738ccaff030SJeremy L Thompson NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 739ccaff030SJeremy L Thompson kilogram = fabs(kilogram); 740ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_Kelvin", 741ccaff030SJeremy L Thompson "1 Kelvin in scaled temperature units", 742ccaff030SJeremy L Thompson NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 743ccaff030SJeremy L Thompson Kelvin = fabs(Kelvin); 744ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 745ccaff030SJeremy L Thompson NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 746ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 747ccaff030SJeremy L Thompson NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 748ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 749ccaff030SJeremy L Thompson NULL, P0, &P0, NULL); CHKERRQ(ierr); 750ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 751ccaff030SJeremy L Thompson NULL, N, &N, NULL); CHKERRQ(ierr); 752ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 753ccaff030SJeremy L Thompson NULL, cv, &cv, NULL); CHKERRQ(ierr); 754ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 755ccaff030SJeremy L Thompson NULL, cp, &cp, NULL); CHKERRQ(ierr); 756ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 757ccaff030SJeremy L Thompson NULL, g, &g, NULL); CHKERRQ(ierr); 758ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lambda", 759ccaff030SJeremy L Thompson "Stokes hypothesis second viscosity coefficient", 760ccaff030SJeremy L Thompson NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 761ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 762ccaff030SJeremy L Thompson NULL, mu, &mu, NULL); CHKERRQ(ierr); 763ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-k", "Thermal conductivity", 764ccaff030SJeremy L Thompson NULL, k, &k, NULL); CHKERRQ(ierr); 765ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-CtauS", 766ccaff030SJeremy L Thompson "Scale coefficient for tau (nondimensional)", 767ccaff030SJeremy L Thompson NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 768ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-strong_form", 769ccaff030SJeremy L Thompson "Strong (1) or weak/integrated by parts (0) advection residual", 770ccaff030SJeremy L Thompson NULL, strong_form, &strong_form, NULL); 771ccaff030SJeremy L Thompson CHKERRQ(ierr); 772ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 773ccaff030SJeremy L Thompson NULL, lx, &lx, NULL); CHKERRQ(ierr); 774ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 775ccaff030SJeremy L Thompson NULL, ly, &ly, NULL); CHKERRQ(ierr); 776ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 777ccaff030SJeremy L Thompson NULL, lz, &lz, NULL); CHKERRQ(ierr); 778ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 779ccaff030SJeremy L Thompson NULL, rc, &rc, NULL); CHKERRQ(ierr); 780ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resx","Target resolution in x", 781ccaff030SJeremy L Thompson NULL, resx, &resx, NULL); CHKERRQ(ierr); 782ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resy","Target resolution in y", 783ccaff030SJeremy L Thompson NULL, resy, &resy, NULL); CHKERRQ(ierr); 784ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resz","Target resolution in z", 785ccaff030SJeremy L Thompson NULL, resz, &resz, NULL); CHKERRQ(ierr); 786ccaff030SJeremy L Thompson PetscInt n = problem->dim; 787ccaff030SJeremy L Thompson ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction", 788ccaff030SJeremy L Thompson NULL, DMBoundaryTypes, (PetscEnum *)periodicity, 789ccaff030SJeremy L Thompson &n, NULL); CHKERRQ(ierr); 790ccaff030SJeremy L Thompson n = problem->dim; 791ccaff030SJeremy L Thompson center[0] = 0.5 * lx; 792ccaff030SJeremy L Thompson center[1] = 0.5 * ly; 793ccaff030SJeremy L Thompson center[2] = 0.5 * lz; 794ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-center", "Location of bubble center", 795ccaff030SJeremy L Thompson NULL, center, &n, NULL); CHKERRQ(ierr); 796ccaff030SJeremy L Thompson n = problem->dim; 797ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-dc_axis", 798ccaff030SJeremy L Thompson "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 799ccaff030SJeremy L Thompson NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 800ccaff030SJeremy L Thompson { 801ccaff030SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 802ccaff030SJeremy L Thompson PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 803ccaff030SJeremy L Thompson if (norm > 0) { 804ccaff030SJeremy L Thompson for (int i=0; i<3; i++) dc_axis[i] /= norm; 805ccaff030SJeremy L Thompson } 806ccaff030SJeremy L Thompson } 807ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-output_freq", 808ccaff030SJeremy L Thompson "Frequency of output, in number of steps", 809ccaff030SJeremy L Thompson NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 810ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-continue", "Continue from previous solution", 811ccaff030SJeremy L Thompson NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 812ff6701fcSJed Brown ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 813ff6701fcSJed Brown NULL, degree, °ree, NULL); CHKERRQ(ierr); 814ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 815ccaff030SJeremy L Thompson NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 816ccaff030SJeremy L Thompson PetscStrncpy(user->outputfolder, ".", 2); 817ccaff030SJeremy L Thompson ierr = PetscOptionsString("-of", "Output folder", 818ccaff030SJeremy L Thompson NULL, user->outputfolder, user->outputfolder, 819ccaff030SJeremy L Thompson sizeof(user->outputfolder), NULL); CHKERRQ(ierr); 820ccaff030SJeremy L Thompson ierr = PetscOptionsEnd(); CHKERRQ(ierr); 821ccaff030SJeremy L Thompson 822ccaff030SJeremy L Thompson // Define derived units 823ccaff030SJeremy L Thompson Pascal = kilogram / (meter * PetscSqr(second)); 824ccaff030SJeremy L Thompson JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 825ccaff030SJeremy L Thompson mpersquareds = meter / PetscSqr(second); 826ccaff030SJeremy L Thompson WpermK = kilogram * meter / (pow(second,3) * Kelvin); 827ccaff030SJeremy L Thompson kgpercubicm = kilogram / pow(meter,3); 828ccaff030SJeremy L Thompson kgpersquaredms = kilogram / (PetscSqr(meter) * second); 829ccaff030SJeremy L Thompson Joulepercubicm = kilogram / (meter * PetscSqr(second)); 830ccaff030SJeremy L Thompson 831ccaff030SJeremy L Thompson // Scale variables to desired units 832ccaff030SJeremy L Thompson theta0 *= Kelvin; 833ccaff030SJeremy L Thompson thetaC *= Kelvin; 834ccaff030SJeremy L Thompson P0 *= Pascal; 835ccaff030SJeremy L Thompson N *= (1./second); 836ccaff030SJeremy L Thompson cv *= JperkgK; 837ccaff030SJeremy L Thompson cp *= JperkgK; 838ccaff030SJeremy L Thompson Rd = cp - cv; 839ccaff030SJeremy L Thompson g *= mpersquareds; 840ccaff030SJeremy L Thompson mu *= Pascal * second; 841ccaff030SJeremy L Thompson k *= WpermK; 842ccaff030SJeremy L Thompson lx = fabs(lx) * meter; 843ccaff030SJeremy L Thompson ly = fabs(ly) * meter; 844ccaff030SJeremy L Thompson lz = fabs(lz) * meter; 845ccaff030SJeremy L Thompson rc = fabs(rc) * meter; 846ccaff030SJeremy L Thompson resx = fabs(resx) * meter; 847ccaff030SJeremy L Thompson resy = fabs(resy) * meter; 848ccaff030SJeremy L Thompson resz = fabs(resz) * meter; 849ccaff030SJeremy L Thompson for (int i=0; i<3; i++) center[i] *= meter; 850ccaff030SJeremy L Thompson 851ccaff030SJeremy L Thompson const CeedInt dim = problem->dim, ncompx = problem->dim, 852ccaff030SJeremy L Thompson qdatasize = problem->qdatasize; 853ccaff030SJeremy L Thompson // Set up the libCEED context 854ccaff030SJeremy L Thompson struct SetupContext_ ctxSetup = { 855ccaff030SJeremy L Thompson .theta0 = theta0, 856ccaff030SJeremy L Thompson .thetaC = thetaC, 857ccaff030SJeremy L Thompson .P0 = P0, 858ccaff030SJeremy L Thompson .N = N, 859ccaff030SJeremy L Thompson .cv = cv, 860ccaff030SJeremy L Thompson .cp = cp, 861ccaff030SJeremy L Thompson .Rd = Rd, 862ccaff030SJeremy L Thompson .g = g, 863ccaff030SJeremy L Thompson .rc = rc, 864ccaff030SJeremy L Thompson .lx = lx, 865ccaff030SJeremy L Thompson .ly = ly, 866ccaff030SJeremy L Thompson .lz = lz, 867ccaff030SJeremy L Thompson .periodicity0 = periodicity[0], 868ccaff030SJeremy L Thompson .periodicity1 = periodicity[1], 869ccaff030SJeremy L Thompson .periodicity2 = periodicity[2], 870ccaff030SJeremy L Thompson .center[0] = center[0], 871ccaff030SJeremy L Thompson .center[1] = center[1], 872ccaff030SJeremy L Thompson .center[2] = center[2], 873ccaff030SJeremy L Thompson .dc_axis[0] = dc_axis[0], 874ccaff030SJeremy L Thompson .dc_axis[1] = dc_axis[1], 875ccaff030SJeremy L Thompson .dc_axis[2] = dc_axis[2], 876ccaff030SJeremy L Thompson .time = 0, 877ccaff030SJeremy L Thompson }; 878ccaff030SJeremy L Thompson 879380e0a58Svaleriabarra // Create the mesh 880ccaff030SJeremy L Thompson { 881ccaff030SJeremy L Thompson const PetscReal scale[3] = {lx, ly, lz}; 882ccaff030SJeremy L Thompson ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 883ccaff030SJeremy L Thompson periodicity, PETSC_TRUE, &dm); 884ccaff030SJeremy L Thompson CHKERRQ(ierr); 885ccaff030SJeremy L Thompson } 886380e0a58Svaleriabarra 887380e0a58Svaleriabarra // Distribute the mesh over processes 888380e0a58Svaleriabarra { 889ccaff030SJeremy L Thompson DM dmDist = NULL; 890ccaff030SJeremy L Thompson PetscPartitioner part; 891ccaff030SJeremy L Thompson 892ccaff030SJeremy L Thompson ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 893ccaff030SJeremy L Thompson ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 894ccaff030SJeremy L Thompson ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 895ccaff030SJeremy L Thompson if (dmDist) { 896ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 897ccaff030SJeremy L Thompson dm = dmDist; 898ccaff030SJeremy L Thompson } 899ccaff030SJeremy L Thompson } 900ccaff030SJeremy L Thompson ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 901ccaff030SJeremy L Thompson 902380e0a58Svaleriabarra // Setup DM 903ccaff030SJeremy L Thompson ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 904ccaff030SJeremy L Thompson ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 905ff6701fcSJed Brown ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr); 906380e0a58Svaleriabarra 907380e0a58Svaleriabarra // Print FEM space information 908ccaff030SJeremy L Thompson if (!test) { 909ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 910380e0a58Svaleriabarra "Degree of FEM space: %D\n", 911ff6701fcSJed Brown degree); CHKERRQ(ierr); 912ccaff030SJeremy L Thompson } 913380e0a58Svaleriabarra 914380e0a58Svaleriabarra // Refine DM for high-order viz 915ccaff030SJeremy L Thompson dmviz = NULL; 916ccaff030SJeremy L Thompson interpviz = NULL; 917ccaff030SJeremy L Thompson if (viz_refine) { 918ff6701fcSJed Brown DM dmhierarchy[viz_refine+1]; 919ff6701fcSJed Brown 920ccaff030SJeremy L Thompson ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 921ff6701fcSJed Brown dmhierarchy[0] = dm; 922ff6701fcSJed Brown for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 923ff6701fcSJed Brown Mat interp_next; 924ff6701fcSJed Brown 925ff6701fcSJed Brown ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 926ccaff030SJeremy L Thompson CHKERRQ(ierr); 927ff6701fcSJed Brown ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 928ff6701fcSJed Brown d = (d + 1) / 2; 929ff6701fcSJed Brown if (i + 1 == viz_refine) d = 1; 930ff6701fcSJed Brown ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr); 931ff6701fcSJed Brown ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 932ff6701fcSJed Brown &interp_next, NULL); CHKERRQ(ierr); 933ff6701fcSJed Brown if (!i) interpviz = interp_next; 934ff6701fcSJed Brown else { 935ff6701fcSJed Brown Mat C; 936ff6701fcSJed Brown ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 937ff6701fcSJed Brown PETSC_DECIDE, &C); CHKERRQ(ierr); 938ff6701fcSJed Brown ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 939ff6701fcSJed Brown ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 940ff6701fcSJed Brown interpviz = C; 941ff6701fcSJed Brown } 942ff6701fcSJed Brown } 943cb3e2689Svaleriabarra for (PetscInt i=1; i<viz_refine; i++) { 944ff6701fcSJed Brown ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 945cb3e2689Svaleriabarra } 946ff6701fcSJed Brown dmviz = dmhierarchy[viz_refine]; 947ccaff030SJeremy L Thompson } 948ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 949ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 950ccaff030SJeremy L Thompson ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 951ccaff030SJeremy L Thompson lnodes /= ncompq; 952ccaff030SJeremy L Thompson 953ccaff030SJeremy L Thompson { 954ccaff030SJeremy L Thompson // Print grid information 955ccaff030SJeremy L Thompson CeedInt gdofs, odofs; 956ccaff030SJeremy L Thompson int comm_size; 957ccaff030SJeremy L Thompson char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 958ccaff030SJeremy L Thompson ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 959ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 960ccaff030SJeremy L Thompson ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 961ccaff030SJeremy L Thompson ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 962ccaff030SJeremy L Thompson sizeof(box_faces_str), NULL); CHKERRQ(ierr); 963ccaff030SJeremy L Thompson if (!test) { 964380e0a58Svaleriabarra ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d rank(s)\n", 965ccaff030SJeremy L Thompson gdofs, odofs, comm_size); CHKERRQ(ierr); 966ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr); 967ccaff030SJeremy L Thompson ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str); 968ccaff030SJeremy L Thompson CHKERRQ(ierr); 969ccaff030SJeremy L Thompson } 970ccaff030SJeremy L Thompson 971ccaff030SJeremy L Thompson } 972ccaff030SJeremy L Thompson 973ccaff030SJeremy L Thompson // Set up global mass vector 974ccaff030SJeremy L Thompson ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 975ccaff030SJeremy L Thompson 976ccaff030SJeremy L Thompson // Set up CEED 977ccaff030SJeremy L Thompson // CEED Bases 978ccaff030SJeremy L Thompson CeedInit(ceedresource, &ceed); 979ccaff030SJeremy L Thompson numP = degree + 1; 980ccaff030SJeremy L Thompson numQ = numP + qextra; 981ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 982ccaff030SJeremy L Thompson &basisq); 983ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 984ccaff030SJeremy L Thompson &basisx); 985ccaff030SJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 986ccaff030SJeremy L Thompson CEED_GAUSS_LOBATTO, &basisxc); 987ccaff030SJeremy L Thompson 988ccaff030SJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 989ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 990ccaff030SJeremy L Thompson CHKERRQ(ierr); 991ccaff030SJeremy L Thompson 992ccaff030SJeremy L Thompson // CEED Restrictions 993ccaff030SJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, degree+1, &restrictq); 994ccaff030SJeremy L Thompson CHKERRQ(ierr); 995ccaff030SJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, &restrictx); CHKERRQ(ierr); 996ccaff030SJeremy L Thompson DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 997ccaff030SJeremy L Thompson localNelem = cEnd - cStart; 998ccaff030SJeremy L Thompson CeedInt numQdim = CeedIntPow(numQ, dim); 999ccaff030SJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, localNelem, numQdim, 1000ccaff030SJeremy L Thompson localNelem*numQdim, qdatasize, 1001ccaff030SJeremy L Thompson CEED_STRIDES_BACKEND, &restrictqdi); 1002ccaff030SJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, localNelem, PetscPowInt(numP, dim), 1003ccaff030SJeremy L Thompson localNelem*PetscPowInt(numP, dim), ncompx, 1004ccaff030SJeremy L Thompson CEED_STRIDES_BACKEND, &restrictxcoord); 1005ccaff030SJeremy L Thompson 1006ccaff030SJeremy L Thompson ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1007ccaff030SJeremy L Thompson ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1008ccaff030SJeremy L Thompson 1009ccaff030SJeremy L Thompson // Create the CEED vectors that will be needed in setup 1010ccaff030SJeremy L Thompson CeedInt Nqpts; 1011ccaff030SJeremy L Thompson CeedBasisGetNumQuadraturePoints(basisq, &Nqpts); 1012ccaff030SJeremy L Thompson CeedVectorCreate(ceed, qdatasize*localNelem*Nqpts, &qdata); 1013ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1014ccaff030SJeremy L Thompson 1015ccaff030SJeremy L Thompson // Create the Q-function that builds the quadrature data for the NS operator 1016ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->setup, problem->setup_loc, 1017ccaff030SJeremy L Thompson &qf_setup); 1018ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_setup, "dx", ncompx*dim, CEED_EVAL_GRAD); 1019ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 1020ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_setup, "qdata", qdatasize, CEED_EVAL_NONE); 1021ccaff030SJeremy L Thompson 1022ccaff030SJeremy L Thompson // Create the Q-function that sets the ICs of the operator 1023ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1024ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1025ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1026ccaff030SJeremy L Thompson 1027ccaff030SJeremy L Thompson qf_rhs = NULL; 1028ccaff030SJeremy L Thompson if (problem->apply_rhs) { // Create the Q-function that defines the action of the RHS operator 1029ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->apply_rhs, 1030ccaff030SJeremy L Thompson problem->apply_rhs_loc, &qf_rhs); 1031ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "q", ncompq, CEED_EVAL_INTERP); 1032ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "dq", ncompq*dim, CEED_EVAL_GRAD); 1033ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "qdata", qdatasize, CEED_EVAL_NONE); 1034ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_rhs, "x", ncompx, CEED_EVAL_INTERP); 1035ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_rhs, "v", ncompq, CEED_EVAL_INTERP); 1036ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_rhs, "dv", ncompq*dim, CEED_EVAL_GRAD); 1037ccaff030SJeremy L Thompson } 1038ccaff030SJeremy L Thompson 1039ccaff030SJeremy L Thompson qf_ifunction = NULL; 1040ccaff030SJeremy L Thompson if (problem->apply_ifunction) { // Create the Q-function that defines the action of the IFunction 1041ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->apply_ifunction, 1042ccaff030SJeremy L Thompson problem->apply_ifunction_loc, &qf_ifunction); 1043ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "q", ncompq, CEED_EVAL_INTERP); 1044ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "dq", ncompq*dim, CEED_EVAL_GRAD); 1045ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "qdot", ncompq, CEED_EVAL_INTERP); 1046ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "qdata", qdatasize, CEED_EVAL_NONE); 1047ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ifunction, "x", ncompx, CEED_EVAL_INTERP); 1048ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ifunction, "v", ncompq, CEED_EVAL_INTERP); 1049ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ifunction, "dv", ncompq*dim, CEED_EVAL_GRAD); 1050ccaff030SJeremy L Thompson } 1051ccaff030SJeremy L Thompson 1052ccaff030SJeremy L Thompson // Create the operator that builds the quadrature data for the NS operator 1053ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup); 1054ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1055ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, 1056ccaff030SJeremy L Thompson basisx, CEED_VECTOR_NONE); 1057ccaff030SJeremy L Thompson CeedOperatorSetField(op_setup, "qdata", restrictqdi, 1058ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1059ccaff030SJeremy L Thompson 1060ccaff030SJeremy L Thompson // Create the operator that sets the ICs 1061ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 1062ccaff030SJeremy L Thompson CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 1063ccaff030SJeremy L Thompson CeedOperatorSetField(op_ics, "q0", restrictq, 1064ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1065ccaff030SJeremy L Thompson 1066ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 1067ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 1068ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1069ccaff030SJeremy L Thompson 1070ccaff030SJeremy L Thompson if (qf_rhs) { // Create the RHS physics operator 1071ccaff030SJeremy L Thompson CeedOperator op; 1072ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op); 1073ccaff030SJeremy L Thompson CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1074ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1075ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdata", restrictqdi, 1076ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 1077ccaff030SJeremy L Thompson CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1078ccaff030SJeremy L Thompson CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1079ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1080ccaff030SJeremy L Thompson user->op_rhs = op; 1081ccaff030SJeremy L Thompson } 1082ccaff030SJeremy L Thompson 1083ccaff030SJeremy L Thompson if (qf_ifunction) { // Create the IFunction operator 1084ccaff030SJeremy L Thompson CeedOperator op; 1085ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ifunction, NULL, NULL, &op); 1086ccaff030SJeremy L Thompson CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1087ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1088ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 1089ccaff030SJeremy L Thompson CeedOperatorSetField(op, "qdata", restrictqdi, 1090ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 1091ccaff030SJeremy L Thompson CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1092ccaff030SJeremy L Thompson CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1093ccaff030SJeremy L Thompson CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1094ccaff030SJeremy L Thompson user->op_ifunction = op; 1095ccaff030SJeremy L Thompson } 1096ccaff030SJeremy L Thompson 1097ccaff030SJeremy L Thompson CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup); 1098ccaff030SJeremy L Thompson CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd}; 1099ccaff030SJeremy L Thompson struct Advection2dContext_ ctxAdvection2d = { 1100ccaff030SJeremy L Thompson .CtauS = CtauS, 1101ccaff030SJeremy L Thompson .strong_form = strong_form, 1102ccaff030SJeremy L Thompson .stabilization = stab, 1103ccaff030SJeremy L Thompson }; 1104ccaff030SJeremy L Thompson switch (problemChoice) { 1105ccaff030SJeremy L Thompson case NS_DENSITY_CURRENT: 1106ccaff030SJeremy L Thompson if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxNS, sizeof ctxNS); 1107ccaff030SJeremy L Thompson if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxNS, 1108ccaff030SJeremy L Thompson sizeof ctxNS); 1109ccaff030SJeremy L Thompson break; 1110ccaff030SJeremy L Thompson case NS_ADVECTION: 1111ccaff030SJeremy L Thompson case NS_ADVECTION2D: 1112ccaff030SJeremy L Thompson if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxAdvection2d, 1113ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 1114ccaff030SJeremy L Thompson if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxAdvection2d, 1115ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 1116ccaff030SJeremy L Thompson } 1117ccaff030SJeremy L Thompson 1118ccaff030SJeremy L Thompson // Set up PETSc context 1119ccaff030SJeremy L Thompson // Set up units structure 1120ccaff030SJeremy L Thompson units->meter = meter; 1121ccaff030SJeremy L Thompson units->kilogram = kilogram; 1122ccaff030SJeremy L Thompson units->second = second; 1123ccaff030SJeremy L Thompson units->Kelvin = Kelvin; 1124ccaff030SJeremy L Thompson units->Pascal = Pascal; 1125ccaff030SJeremy L Thompson units->JperkgK = JperkgK; 1126ccaff030SJeremy L Thompson units->mpersquareds = mpersquareds; 1127ccaff030SJeremy L Thompson units->WpermK = WpermK; 1128ccaff030SJeremy L Thompson units->kgpercubicm = kgpercubicm; 1129ccaff030SJeremy L Thompson units->kgpersquaredms = kgpersquaredms; 1130ccaff030SJeremy L Thompson units->Joulepercubicm = Joulepercubicm; 1131ccaff030SJeremy L Thompson 1132ccaff030SJeremy L Thompson // Set up user structure 1133ccaff030SJeremy L Thompson user->comm = comm; 1134ccaff030SJeremy L Thompson user->outputfreq = outputfreq; 1135ccaff030SJeremy L Thompson user->contsteps = contsteps; 1136ccaff030SJeremy L Thompson user->units = units; 1137ccaff030SJeremy L Thompson user->dm = dm; 1138ccaff030SJeremy L Thompson user->dmviz = dmviz; 1139ccaff030SJeremy L Thompson user->interpviz = interpviz; 1140ccaff030SJeremy L Thompson user->ceed = ceed; 1141ccaff030SJeremy L Thompson 1142ccaff030SJeremy L Thompson // Calculate qdata and ICs 1143ccaff030SJeremy L Thompson // Set up state global and local vectors 1144ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1145ccaff030SJeremy L Thompson 1146ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1147ccaff030SJeremy L Thompson 1148ccaff030SJeremy L Thompson // Apply Setup Ceed Operators 1149ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 1150ccaff030SJeremy L Thompson CeedOperatorApply(op_setup, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 1151ccaff030SJeremy L Thompson ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1152ccaff030SJeremy L Thompson user->M); CHKERRQ(ierr); 1153ccaff030SJeremy L Thompson 115427dd567eSvaleriabarra ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 1155ccaff030SJeremy L Thompson &ctxSetup, 0.0); 1156ccaff030SJeremy L Thompson if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1157ccaff030SJeremy L Thompson // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1158ccaff030SJeremy L Thompson // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1159ccaff030SJeremy L Thompson // slower execution. 1160ccaff030SJeremy L Thompson Vec Qbc; 1161ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1162ccaff030SJeremy L Thompson ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1163ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1164ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1165ccaff030SJeremy L Thompson ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1166ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1167ccaff030SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)dm, 1168380e0a58Svaleriabarra "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 1169380e0a58Svaleriabarra CHKERRQ(ierr); 1170ccaff030SJeremy L Thompson } 1171ccaff030SJeremy L Thompson 1172ccaff030SJeremy L Thompson MPI_Comm_rank(comm, &rank); 1173ccaff030SJeremy L Thompson if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);} 1174ccaff030SJeremy L Thompson // Gather initial Q values 1175ccaff030SJeremy L Thompson // In case of continuation of simulation, set up initial values from binary file 1176ccaff030SJeremy L Thompson if (contsteps) { // continue from existent solution 1177ccaff030SJeremy L Thompson PetscViewer viewer; 1178ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1179ccaff030SJeremy L Thompson // Read input 1180ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1181ccaff030SJeremy L Thompson user->outputfolder); 1182ccaff030SJeremy L Thompson CHKERRQ(ierr); 1183ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1184ccaff030SJeremy L Thompson CHKERRQ(ierr); 1185ccaff030SJeremy L Thompson ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1186ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1187ccaff030SJeremy L Thompson } 1188ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1189ccaff030SJeremy L Thompson 1190ccaff030SJeremy L Thompson // Create and setup TS 1191ccaff030SJeremy L Thompson ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1192ccaff030SJeremy L Thompson ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1193ccaff030SJeremy L Thompson if (implicit) { 1194ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1195ccaff030SJeremy L Thompson if (user->op_ifunction) { 1196ccaff030SJeremy L Thompson ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1197ccaff030SJeremy L Thompson } else { // Implicit integrators can fall back to using an RHSFunction 1198ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1199ccaff030SJeremy L Thompson } 1200ccaff030SJeremy L Thompson } else { 1201ccaff030SJeremy L Thompson if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL, 1202ccaff030SJeremy L Thompson "Problem does not provide RHSFunction"); 1203ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1204ccaff030SJeremy L Thompson ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1205ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1206ccaff030SJeremy L Thompson } 1207ccaff030SJeremy L Thompson ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1208ccaff030SJeremy L Thompson ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1209ccaff030SJeremy L Thompson ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 1210ccaff030SJeremy L Thompson if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);} 1211ccaff030SJeremy L Thompson ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1212ccaff030SJeremy L Thompson ierr = TSAdaptSetStepLimits(adapt, 1213ccaff030SJeremy L Thompson 1.e-12 * units->second, 1214ccaff030SJeremy L Thompson 1.e2 * units->second); CHKERRQ(ierr); 1215ccaff030SJeremy L Thompson ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1216ccaff030SJeremy L Thompson if (!contsteps) { // print initial condition 1217ccaff030SJeremy L Thompson if (!test) { 1218ccaff030SJeremy L Thompson ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1219ccaff030SJeremy L Thompson } 1220ccaff030SJeremy L Thompson } else { // continue from time of last output 1221ccaff030SJeremy L Thompson PetscReal time; 1222ccaff030SJeremy L Thompson PetscInt count; 1223ccaff030SJeremy L Thompson PetscViewer viewer; 1224ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1225ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1226ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 1227ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1228ccaff030SJeremy L Thompson CHKERRQ(ierr); 1229ccaff030SJeremy L Thompson ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1230ccaff030SJeremy L Thompson CHKERRQ(ierr); 1231ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1232ccaff030SJeremy L Thompson ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1233ccaff030SJeremy L Thompson } 1234ccaff030SJeremy L Thompson if (!test) { 1235ccaff030SJeremy L Thompson ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1236ccaff030SJeremy L Thompson } 1237ccaff030SJeremy L Thompson 1238ccaff030SJeremy L Thompson // Solve 1239ccaff030SJeremy L Thompson start = MPI_Wtime(); 1240ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1241ccaff030SJeremy L Thompson ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1242ccaff030SJeremy L Thompson cpu_time_used = MPI_Wtime() - start; 1243ccaff030SJeremy L Thompson ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1244ccaff030SJeremy L Thompson ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1245ccaff030SJeremy L Thompson comm); CHKERRQ(ierr); 1246ccaff030SJeremy L Thompson if (!test) { 1247ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1248380e0a58Svaleriabarra "Time taken for solution (sec): %g\n", 1249ccaff030SJeremy L Thompson (double)cpu_time_used); CHKERRQ(ierr); 1250ccaff030SJeremy L Thompson } 1251ccaff030SJeremy L Thompson 1252ccaff030SJeremy L Thompson // Get error 1253ccaff030SJeremy L Thompson if (problem->non_zero_time && !test) { 1254ccaff030SJeremy L Thompson Vec Qexact, Qexactloc; 1255ccaff030SJeremy L Thompson PetscReal norm; 1256ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1257ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1258ccaff030SJeremy L Thompson ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1259ccaff030SJeremy L Thompson 126027dd567eSvaleriabarra ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 1261ccaff030SJeremy L Thompson restrictq, &ctxSetup, ftime); CHKERRQ(ierr); 1262ccaff030SJeremy L Thompson 1263ccaff030SJeremy L Thompson ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1264ccaff030SJeremy L Thompson ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1265ccaff030SJeremy L Thompson CeedVectorDestroy(&q0ceed); 1266ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1267ccaff030SJeremy L Thompson "Max Error: %g\n", 1268ccaff030SJeremy L Thompson (double)norm); CHKERRQ(ierr); 1269357222b3Svaleriabarra // Clean up vectors 1270357222b3Svaleriabarra ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1271357222b3Svaleriabarra ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1272ccaff030SJeremy L Thompson } 1273ccaff030SJeremy L Thompson 1274ccaff030SJeremy L Thompson // Output Statistics 1275ccaff030SJeremy L Thompson ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr); 1276ccaff030SJeremy L Thompson if (!test) { 1277ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1278ccaff030SJeremy L Thompson "Time integrator took %D time steps to reach final time %g\n", 1279ccaff030SJeremy L Thompson steps, (double)ftime); CHKERRQ(ierr); 1280ccaff030SJeremy L Thompson } 1281ccaff030SJeremy L Thompson 1282ccaff030SJeremy L Thompson // Clean up libCEED 1283ccaff030SJeremy L Thompson CeedVectorDestroy(&qdata); 1284ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qceed); 1285ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qdotceed); 1286ccaff030SJeremy L Thompson CeedVectorDestroy(&user->gceed); 1287ccaff030SJeremy L Thompson CeedVectorDestroy(&xcorners); 1288ccaff030SJeremy L Thompson CeedBasisDestroy(&basisq); 1289ccaff030SJeremy L Thompson CeedBasisDestroy(&basisx); 1290ccaff030SJeremy L Thompson CeedBasisDestroy(&basisxc); 1291ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictq); 1292ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictx); 1293ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictqdi); 1294ccaff030SJeremy L Thompson CeedElemRestrictionDestroy(&restrictxcoord); 1295ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_setup); 1296ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ics); 1297ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_rhs); 1298ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ifunction); 1299ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_setup); 1300ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_ics); 1301ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_rhs); 1302ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_ifunction); 1303ccaff030SJeremy L Thompson CeedDestroy(&ceed); 1304ccaff030SJeremy L Thompson 1305ccaff030SJeremy L Thompson // Clean up PETSc 1306ccaff030SJeremy L Thompson ierr = VecDestroy(&Q); CHKERRQ(ierr); 1307ccaff030SJeremy L Thompson ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1308ccaff030SJeremy L Thompson ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1309ccaff030SJeremy L Thompson ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1310ccaff030SJeremy L Thompson ierr = TSDestroy(&ts); CHKERRQ(ierr); 1311ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1312ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 1313ccaff030SJeremy L Thompson ierr = PetscFree(user); CHKERRQ(ierr); 1314ccaff030SJeremy L Thompson return PetscFinalize(); 1315ccaff030SJeremy L Thompson } 1316ccaff030SJeremy L Thompson 1317