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 // 3184d34d69SLeila Ghaffari // ./navierstokes -ceed /cpu/self -problem density_current -degree 1 3284d34d69SLeila Ghaffari // ./navierstokes -ceed /gpu/occa -problem advection -degree 1 33ccaff030SJeremy L Thompson // 3484d34d69SLeila Ghaffari //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 3584d34d69SLeila Ghaffari //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 3684d34d69SLeila Ghaffari //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" 49b0137797SLeila Ghaffari #include "setup-boundary.h" 50ccaff030SJeremy L Thompson #include "advection.h" 51ccaff030SJeremy L Thompson #include "advection2d.h" 52ccaff030SJeremy L Thompson #include "densitycurrent.h" 53ccaff030SJeremy L Thompson 5484d34d69SLeila Ghaffari #if PETSC_VERSION_LT(3,14,0) 5584d34d69SLeila Ghaffari # define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i) 5684d34d69SLeila Ghaffari # define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i) 5784d34d69SLeila Ghaffari #endif 5884d34d69SLeila Ghaffari 5984d34d69SLeila Ghaffari // MemType Options 6084d34d69SLeila Ghaffari static const char *const memTypes[] = { 6184d34d69SLeila Ghaffari "host", 6284d34d69SLeila Ghaffari "device", 6384d34d69SLeila Ghaffari "memType", "CEED_MEM_", NULL 6484d34d69SLeila Ghaffari }; 6584d34d69SLeila Ghaffari 66ccaff030SJeremy L Thompson // Problem Options 67ccaff030SJeremy L Thompson typedef enum { 68ccaff030SJeremy L Thompson NS_DENSITY_CURRENT = 0, 69ccaff030SJeremy L Thompson NS_ADVECTION = 1, 70ccaff030SJeremy L Thompson NS_ADVECTION2D = 2, 71ccaff030SJeremy L Thompson } problemType; 72ccaff030SJeremy L Thompson static const char *const problemTypes[] = { 73ccaff030SJeremy L Thompson "density_current", 74ccaff030SJeremy L Thompson "advection", 75ccaff030SJeremy L Thompson "advection2d", 7684d34d69SLeila Ghaffari "problemType", "NS_", NULL 77ccaff030SJeremy L Thompson }; 78ccaff030SJeremy L Thompson 791184866aSLeila Ghaffari // Wind Options for Advection 801184866aSLeila Ghaffari typedef enum { 811184866aSLeila Ghaffari ADVECTION_WIND_ROTATION = 0, 821184866aSLeila Ghaffari ADVECTION_WIND_TRANSLATION = 1, 831184866aSLeila Ghaffari } WindType; 841184866aSLeila Ghaffari static const char *const WindTypes[] = { 851184866aSLeila Ghaffari "rotation", 861184866aSLeila Ghaffari "translation", 871184866aSLeila Ghaffari "WindType", "ADVECTION_WIND_", NULL 881184866aSLeila Ghaffari }; 891184866aSLeila Ghaffari 90ccaff030SJeremy L Thompson typedef enum { 91ccaff030SJeremy L Thompson STAB_NONE = 0, 92ccaff030SJeremy L Thompson STAB_SU = 1, // Streamline Upwind 93ccaff030SJeremy L Thompson STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 94ccaff030SJeremy L Thompson } StabilizationType; 95ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = { 9684d34d69SLeila Ghaffari "none", 97ccaff030SJeremy L Thompson "SU", 98ccaff030SJeremy L Thompson "SUPG", 99ccaff030SJeremy L Thompson "StabilizationType", "STAB_", NULL 100ccaff030SJeremy L Thompson }; 101ccaff030SJeremy L Thompson 10284d34d69SLeila Ghaffari // Test Options 10384d34d69SLeila Ghaffari typedef enum { 10484d34d69SLeila Ghaffari TEST_NONE = 0, // Non test mode 10584d34d69SLeila Ghaffari TEST_EXPLICIT = 1, // Explicit test 10684d34d69SLeila Ghaffari TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab 10784d34d69SLeila Ghaffari TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab 10884d34d69SLeila Ghaffari } testType; 10984d34d69SLeila Ghaffari static const char *const testTypes[] = { 11084d34d69SLeila Ghaffari "none", 11184d34d69SLeila Ghaffari "explicit", 11284d34d69SLeila Ghaffari "implicit_stab_none", 11384d34d69SLeila Ghaffari "implicit_stab_supg", 11484d34d69SLeila Ghaffari "testType", "TEST_", NULL 11584d34d69SLeila Ghaffari }; 11684d34d69SLeila Ghaffari 11784d34d69SLeila Ghaffari // Tests specific data 11884d34d69SLeila Ghaffari typedef struct { 11984d34d69SLeila Ghaffari PetscScalar testtol; 12084d34d69SLeila Ghaffari const char *filepath; 12184d34d69SLeila Ghaffari } testData; 12284d34d69SLeila Ghaffari 12384d34d69SLeila Ghaffari testData testOptions[] = { 12484d34d69SLeila Ghaffari [TEST_NONE] = { 12584d34d69SLeila Ghaffari .testtol = 0., 12684d34d69SLeila Ghaffari .filepath = NULL 12784d34d69SLeila Ghaffari }, 12884d34d69SLeila Ghaffari [TEST_EXPLICIT] = { 12984d34d69SLeila Ghaffari .testtol = 1E-5, 13084d34d69SLeila Ghaffari .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin" 13184d34d69SLeila Ghaffari }, 13284d34d69SLeila Ghaffari [TEST_IMPLICIT_STAB_NONE] = { 13384d34d69SLeila Ghaffari .testtol = 5E-4, 13484d34d69SLeila Ghaffari .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin" 13584d34d69SLeila Ghaffari }, 13684d34d69SLeila Ghaffari [TEST_IMPLICIT_STAB_SUPG] = { 13784d34d69SLeila Ghaffari .testtol = 5E-4, 13884d34d69SLeila Ghaffari .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin" 13984d34d69SLeila Ghaffari } 14084d34d69SLeila Ghaffari }; 14184d34d69SLeila Ghaffari 142ccaff030SJeremy L Thompson // Problem specific data 143ccaff030SJeremy L Thompson typedef struct { 1448b982baeSLeila Ghaffari CeedInt dim, qdatasizeVol, qdatasizeSur; 1455603407fSLeila Ghaffari CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction, 1469fe13df9SLeila Ghaffari applyOut_rhs, applyOut_ifunction, applyIn_rhs, applyIn_ifunction; 147ccaff030SJeremy L Thompson PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 148ccaff030SJeremy L Thompson PetscScalar[], void *); 1495603407fSLeila Ghaffari const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc, *applyVol_ifunction_loc, 1509fe13df9SLeila Ghaffari *applyOut_rhs_loc, *applyOut_ifunction_loc, *applyIn_rhs_loc, *applyIn_ifunction_loc; 151ccaff030SJeremy L Thompson const bool non_zero_time; 152ccaff030SJeremy L Thompson } problemData; 153ccaff030SJeremy L Thompson 154ccaff030SJeremy L Thompson problemData problemOptions[] = { 155ccaff030SJeremy L Thompson [NS_DENSITY_CURRENT] = { 156ccaff030SJeremy L Thompson .dim = 3, 157ea6e0f84SLeila Ghaffari .qdatasizeVol = 10, 1588b982baeSLeila Ghaffari .qdatasizeSur = 4, 159b0137797SLeila Ghaffari .setupVol = Setup, 160b0137797SLeila Ghaffari .setupVol_loc = Setup_loc, 161356fbf4bSLeila Ghaffari .setupSur = SetupBoundary, 162356fbf4bSLeila Ghaffari .setupSur_loc = SetupBoundary_loc, 163ccaff030SJeremy L Thompson .ics = ICsDC, 164ccaff030SJeremy L Thompson .ics_loc = ICsDC_loc, 165c96c872fSLeila Ghaffari .applyVol_rhs = DC, 166c96c872fSLeila Ghaffari .applyVol_rhs_loc = DC_loc, 167c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_DC, 168c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_DC_loc, 169ccaff030SJeremy L Thompson .bc = Exact_DC, 17084d34d69SLeila Ghaffari .non_zero_time = PETSC_FALSE, 171ccaff030SJeremy L Thompson }, 172ccaff030SJeremy L Thompson [NS_ADVECTION] = { 173ccaff030SJeremy L Thompson .dim = 3, 174ea6e0f84SLeila Ghaffari .qdatasizeVol = 10, 1758b982baeSLeila Ghaffari .qdatasizeSur = 4, 176b0137797SLeila Ghaffari .setupVol = Setup, 177b0137797SLeila Ghaffari .setupVol_loc = Setup_loc, 178356fbf4bSLeila Ghaffari .setupSur = SetupBoundary, 179356fbf4bSLeila Ghaffari .setupSur_loc = SetupBoundary_loc, 180ccaff030SJeremy L Thompson .ics = ICsAdvection, 181ccaff030SJeremy L Thompson .ics_loc = ICsAdvection_loc, 182c96c872fSLeila Ghaffari .applyVol_rhs = Advection, 183c96c872fSLeila Ghaffari .applyVol_rhs_loc = Advection_loc, 1845603407fSLeila Ghaffari .applyOut_rhs = Advection_Out, 1855603407fSLeila Ghaffari .applyOut_rhs_loc = Advection_Out_loc, 1869fe13df9SLeila Ghaffari .applyIn_rhs = Advection_In, 1879fe13df9SLeila Ghaffari .applyIn_rhs_loc = Advection_In_loc, 188c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_Advection, 189c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_Advection_loc, 1905603407fSLeila Ghaffari .applyOut_ifunction = IFunction_Advection_Out, 1915603407fSLeila Ghaffari .applyOut_ifunction_loc = IFunction_Advection_Out_loc, 1929fe13df9SLeila Ghaffari .applyIn_ifunction = IFunction_Advection_In, 1939fe13df9SLeila Ghaffari .applyIn_ifunction_loc = IFunction_Advection_In_loc, 194ccaff030SJeremy L Thompson .bc = Exact_Advection, 19584d34d69SLeila Ghaffari .non_zero_time = PETSC_FALSE, 196ccaff030SJeremy L Thompson }, 197ccaff030SJeremy L Thompson [NS_ADVECTION2D] = { 198ccaff030SJeremy L Thompson .dim = 2, 199ea6e0f84SLeila Ghaffari .qdatasizeVol = 5, 2008b982baeSLeila Ghaffari .qdatasizeSur = 3, 201c96c872fSLeila Ghaffari .setupVol = Setup2d, 202c96c872fSLeila Ghaffari .setupVol_loc = Setup2d_loc, 203b0137797SLeila Ghaffari .setupSur = SetupBoundary2d, 204b0137797SLeila Ghaffari .setupSur_loc = SetupBoundary2d_loc, 205ccaff030SJeremy L Thompson .ics = ICsAdvection2d, 206ccaff030SJeremy L Thompson .ics_loc = ICsAdvection2d_loc, 207c96c872fSLeila Ghaffari .applyVol_rhs = Advection2d, 208c96c872fSLeila Ghaffari .applyVol_rhs_loc = Advection2d_loc, 2095603407fSLeila Ghaffari .applyOut_rhs = Advection2d_Out, 2105603407fSLeila Ghaffari .applyOut_rhs_loc = Advection2d_Out_loc, 2119fe13df9SLeila Ghaffari .applyIn_rhs = Advection2d_In, 2129fe13df9SLeila Ghaffari .applyIn_rhs_loc = Advection2d_In_loc, 213c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_Advection2d, 214c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_Advection2d_loc, 2155603407fSLeila Ghaffari .applyOut_ifunction = IFunction_Advection2d_Out, 2165603407fSLeila Ghaffari .applyOut_ifunction_loc = IFunction_Advection2d_Out_loc, 2179fe13df9SLeila Ghaffari .applyIn_ifunction = IFunction_Advection2d_In, 2189fe13df9SLeila Ghaffari .applyIn_ifunction_loc = IFunction_Advection2d_In_loc, 219ccaff030SJeremy L Thompson .bc = Exact_Advection2d, 22084d34d69SLeila Ghaffari .non_zero_time = PETSC_TRUE, 221ccaff030SJeremy L Thompson }, 222ccaff030SJeremy L Thompson }; 223ccaff030SJeremy L Thompson 224ccaff030SJeremy L Thompson // PETSc user data 225ccaff030SJeremy L Thompson typedef struct User_ *User; 226ccaff030SJeremy L Thompson typedef struct Units_ *Units; 227ccaff030SJeremy L Thompson 228ccaff030SJeremy L Thompson struct User_ { 229ccaff030SJeremy L Thompson MPI_Comm comm; 230ccaff030SJeremy L Thompson PetscInt outputfreq; 231ccaff030SJeremy L Thompson DM dm; 232ccaff030SJeremy L Thompson DM dmviz; 233ccaff030SJeremy L Thompson Mat interpviz; 234ccaff030SJeremy L Thompson Ceed ceed; 235ccaff030SJeremy L Thompson Units units; 236ccaff030SJeremy L Thompson CeedVector qceed, qdotceed, gceed; 2371e150236SLeila Ghaffari CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction; 238ccaff030SJeremy L Thompson Vec M; 239ccaff030SJeremy L Thompson char outputfolder[PETSC_MAX_PATH_LEN]; 240ccaff030SJeremy L Thompson PetscInt contsteps; 241ccaff030SJeremy L Thompson }; 242ccaff030SJeremy L Thompson 243ccaff030SJeremy L Thompson struct Units_ { 244ccaff030SJeremy L Thompson // fundamental units 245ccaff030SJeremy L Thompson PetscScalar meter; 246ccaff030SJeremy L Thompson PetscScalar kilogram; 247ccaff030SJeremy L Thompson PetscScalar second; 248ccaff030SJeremy L Thompson PetscScalar Kelvin; 249ccaff030SJeremy L Thompson // derived units 250ccaff030SJeremy L Thompson PetscScalar Pascal; 251ccaff030SJeremy L Thompson PetscScalar JperkgK; 252ccaff030SJeremy L Thompson PetscScalar mpersquareds; 253ccaff030SJeremy L Thompson PetscScalar WpermK; 254ccaff030SJeremy L Thompson PetscScalar kgpercubicm; 255ccaff030SJeremy L Thompson PetscScalar kgpersquaredms; 256ccaff030SJeremy L Thompson PetscScalar Joulepercubicm; 257ccaff030SJeremy L Thompson }; 258ccaff030SJeremy L Thompson 259ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC; 260ccaff030SJeremy L Thompson struct SimpleBC_ { 2619fe13df9SLeila Ghaffari PetscInt nwall, nslip[3], noutflow, ninflow; 2629fe13df9SLeila Ghaffari PetscInt walls[6], slips[3][6], outflow[6], inflow[6]; 26384d34d69SLeila Ghaffari PetscBool userbc; 264ccaff030SJeremy L Thompson }; 265ccaff030SJeremy L Thompson 266ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1). 267ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) { 268ccaff030SJeremy L Thompson return i >= 0 ? i : -(i+1); 269ccaff030SJeremy L Thompson } 270ccaff030SJeremy L Thompson 271ccaff030SJeremy L Thompson // Utility function to create local CEED restriction 272ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 27384d34d69SLeila Ghaffari CeedInt height, DMLabel domainLabel, CeedInt value, 274ccaff030SJeremy L Thompson CeedElemRestriction *Erestrict) { 275ccaff030SJeremy L Thompson 276ccaff030SJeremy L Thompson PetscSection section; 2771184866aSLeila Ghaffari PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth; 2780c6c0b13SLeila Ghaffari DMLabel depthLabel; 2790c6c0b13SLeila Ghaffari IS depthIS, iterIS; 28084d34d69SLeila Ghaffari Vec Uloc; 2810c6c0b13SLeila Ghaffari const PetscInt *iterIndices; 282ccaff030SJeremy L Thompson PetscErrorCode ierr; 283ccaff030SJeremy L Thompson 284ccaff030SJeremy L Thompson PetscFunctionBeginUser; 285ccaff030SJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 286da51bdd9SLeila Ghaffari dim -= height; 287ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 288ccaff030SJeremy L Thompson ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 289ccaff030SJeremy L Thompson PetscInt ncomp[nfields], fieldoff[nfields+1]; 290ccaff030SJeremy L Thompson fieldoff[0] = 0; 291ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 292ccaff030SJeremy L Thompson ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 293ccaff030SJeremy L Thompson fieldoff[f+1] = fieldoff[f] + ncomp[f]; 294ccaff030SJeremy L Thompson } 295ccaff030SJeremy L Thompson 2960c6c0b13SLeila Ghaffari ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 2970c6c0b13SLeila Ghaffari ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr); 2980c6c0b13SLeila Ghaffari ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr); 2990c6c0b13SLeila Ghaffari if (domainLabel) { 3000c6c0b13SLeila Ghaffari IS domainIS; 3010c6c0b13SLeila Ghaffari ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr); 3020c6c0b13SLeila Ghaffari ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr); 3030c6c0b13SLeila Ghaffari ierr = ISDestroy(&domainIS); CHKERRQ(ierr); 3040c6c0b13SLeila Ghaffari ierr = ISDestroy(&depthIS); CHKERRQ(ierr); 3050c6c0b13SLeila Ghaffari } else { 3060c6c0b13SLeila Ghaffari iterIS = depthIS; 3070c6c0b13SLeila Ghaffari } 3080c6c0b13SLeila Ghaffari ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr); 3090c6c0b13SLeila Ghaffari ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr); 310ccaff030SJeremy L Thompson ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 3110c6c0b13SLeila Ghaffari for (p=0,eoffset=0; p<Nelem; p++) { 3120c6c0b13SLeila Ghaffari PetscInt c = iterIndices[p]; 313ccaff030SJeremy L Thompson PetscInt numindices, *indices, nnodes; 31484d34d69SLeila Ghaffari ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 31584d34d69SLeila Ghaffari &numindices, &indices, NULL, NULL); 31684d34d69SLeila Ghaffari CHKERRQ(ierr); 31732b5ec5fSJed Brown bool flip = false; 31832b5ec5fSJed Brown if (height > 0) { 31932b5ec5fSJed Brown PetscInt numCells, numFaces, start = -1; 32032b5ec5fSJed Brown const PetscInt *orients, *faces, *cells; 32132b5ec5fSJed Brown ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 32232b5ec5fSJed Brown ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr); 32332b5ec5fSJed Brown if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Expected one cell in support of exterior face, but got %D cells", numCells); 32432b5ec5fSJed Brown ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 32532b5ec5fSJed Brown ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr); 32632b5ec5fSJed Brown for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;} 32732b5ec5fSJed Brown if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Could not find face %D in cone of its support", c); 32832b5ec5fSJed Brown ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 32932b5ec5fSJed Brown if (orients[start] < 0) flip = true; 33032b5ec5fSJed Brown } 33184d34d69SLeila Ghaffari if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 33284d34d69SLeila Ghaffari PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 33384d34d69SLeila Ghaffari c); 334ccaff030SJeremy L Thompson nnodes = numindices / fieldoff[nfields]; 335ccaff030SJeremy L Thompson for (PetscInt i=0; i<nnodes; i++) { 33632b5ec5fSJed Brown PetscInt ii = i; 33732b5ec5fSJed Brown if (flip) { 33832b5ec5fSJed Brown if (P == nnodes) ii = nnodes - 1 - i; 33932b5ec5fSJed Brown else if (P*P == nnodes) { 34032b5ec5fSJed Brown PetscInt row = i / P, col = i % P; 34132b5ec5fSJed Brown ii = row + col * P; 34232b5ec5fSJed Brown } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for flipping point with %D nodes != P (%D) or P^2", nnodes, P); 34332b5ec5fSJed Brown } 344ccaff030SJeremy L Thompson // Check that indices are blocked by node and thus can be coalesced as a single field with 345ccaff030SJeremy L Thompson // fieldoff[nfields] = sum(ncomp) components. 346ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 347ccaff030SJeremy L Thompson for (PetscInt j=0; j<ncomp[f]; j++) { 34832b5ec5fSJed Brown if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j]) 34932b5ec5fSJed Brown != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j) 350ccaff030SJeremy L Thompson SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 351ccaff030SJeremy L Thompson "Cell %D closure indices not interlaced for node %D field %D component %D", 35232b5ec5fSJed Brown c, ii, f, j); 353ccaff030SJeremy L Thompson } 354ccaff030SJeremy L Thompson } 355ccaff030SJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 35632b5ec5fSJed Brown PetscInt loc = Involute(indices[ii*ncomp[0]]); 3576f55dfd5Svaleriabarra erestrict[eoffset++] = loc; 358ccaff030SJeremy L Thompson } 35984d34d69SLeila Ghaffari ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 36084d34d69SLeila Ghaffari &numindices, &indices, NULL, NULL); 36184d34d69SLeila Ghaffari CHKERRQ(ierr); 362ccaff030SJeremy L Thompson } 3630c6c0b13SLeila Ghaffari if (eoffset != Nelem*PetscPowInt(P, dim)) 3640c6c0b13SLeila Ghaffari SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 3650c6c0b13SLeila Ghaffari "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 366ccaff030SJeremy L Thompson PetscPowInt(P, dim),eoffset); 3670c6c0b13SLeila Ghaffari ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr); 3680c6c0b13SLeila Ghaffari ierr = ISDestroy(&iterIS); CHKERRQ(ierr); 3690c6c0b13SLeila Ghaffari 370ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 371ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 372ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 3736f55dfd5Svaleriabarra CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields], 3746f55dfd5Svaleriabarra 1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 3756f55dfd5Svaleriabarra Erestrict); 376ccaff030SJeremy L Thompson ierr = PetscFree(erestrict); CHKERRQ(ierr); 377ccaff030SJeremy L Thompson PetscFunctionReturn(0); 378ccaff030SJeremy L Thompson } 379ccaff030SJeremy L Thompson 380c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain 3811e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, 3821e150236SLeila Ghaffari DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize, 3831e150236SLeila Ghaffari CeedElemRestriction *restrictq, CeedElemRestriction *restrictx, 3841e150236SLeila Ghaffari CeedElemRestriction *restrictqdi) { 385c96c872fSLeila Ghaffari 386c96c872fSLeila Ghaffari DM dmcoord; 3871e150236SLeila Ghaffari CeedInt dim, localNelem; 3881e150236SLeila Ghaffari CeedInt Qdim; 389c96c872fSLeila Ghaffari PetscErrorCode ierr; 390c96c872fSLeila Ghaffari 391c96c872fSLeila Ghaffari PetscFunctionBeginUser; 3921e150236SLeila Ghaffari ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 3931e150236SLeila Ghaffari dim -= height; 3941e150236SLeila Ghaffari Qdim = CeedIntPow(Q, dim); 395c96c872fSLeila Ghaffari ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 396c96c872fSLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 397c96c872fSLeila Ghaffari CHKERRQ(ierr); 398c96c872fSLeila Ghaffari ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq); 399c96c872fSLeila Ghaffari CHKERRQ(ierr); 400c96c872fSLeila Ghaffari ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx); 401c96c872fSLeila Ghaffari CHKERRQ(ierr); 402c96c872fSLeila Ghaffari CeedElemRestrictionGetNumElements(*restrictq, &localNelem); 403c96c872fSLeila Ghaffari CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim, 404c96c872fSLeila Ghaffari qdatasize, qdatasize*localNelem*Qdim, 405c96c872fSLeila Ghaffari CEED_STRIDES_BACKEND, restrictqdi); 406c96c872fSLeila Ghaffari PetscFunctionReturn(0); 407c96c872fSLeila Ghaffari } 408c96c872fSLeila Ghaffari 4091e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain 4101e150236SLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, CeedOperator op_applyVol, 4119fe13df9SLeila Ghaffari CeedQFunction qf_applyOut, CeedQFunction qf_applyIn, CeedQFunction qf_setupSur, CeedInt height, 4129fe13df9SLeila Ghaffari PetscInt nOut, PetscInt valueOut[6], PetscInt nIn, PetscInt valueIn[6], CeedInt numP_Sur, 4139fe13df9SLeila Ghaffari CeedInt numQ_Sur, CeedInt qdatasizeSur, CeedInt NqptsSur, 4141e150236SLeila Ghaffari CeedBasis basisxSur, CeedBasis basisqSur, CeedOperator *op_apply) { 415ca3ac6ddSLeila Ghaffari 4169fe13df9SLeila Ghaffari CeedElemRestriction restrictxOut[6], restrictqOut[6], restrictqdiOut[6], 4179fe13df9SLeila Ghaffari restrictxIn[6], restrictqIn[6], restrictqdiIn[6]; 4181184866aSLeila Ghaffari PetscInt lsize, localNelemOut[6], localNelemIn[6]; 4191e150236SLeila Ghaffari Vec Xloc; 4209fe13df9SLeila Ghaffari CeedVector xcorners, qdataOut[6], qdataIn[6]; 4219fe13df9SLeila Ghaffari CeedOperator op_setupOut[6], op_applyOut[6], op_setupIn[6], op_applyIn[6]; 422ca3ac6ddSLeila Ghaffari DMLabel domainLabel; 4231e150236SLeila Ghaffari PetscScalar *x; 424ca3ac6ddSLeila Ghaffari PetscErrorCode ierr; 425ca3ac6ddSLeila Ghaffari 426ca3ac6ddSLeila Ghaffari PetscFunctionBeginUser; 427ca3ac6ddSLeila Ghaffari // Composite Operaters 428ca3ac6ddSLeila Ghaffari CeedCompositeOperatorCreate(ceed, op_apply); 4291e150236SLeila Ghaffari CeedCompositeOperatorAddSub(*op_apply, op_applyVol); // Apply a Sub-Operator for the volume 430ca3ac6ddSLeila Ghaffari 4319fe13df9SLeila Ghaffari if (nOut || nIn) { 4321e150236SLeila Ghaffari ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 4331e150236SLeila Ghaffari ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr); 4341e150236SLeila Ghaffari ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr); 4351e150236SLeila Ghaffari ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr); 4361e150236SLeila Ghaffari CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x); 437ca3ac6ddSLeila Ghaffari ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr); 4389fe13df9SLeila Ghaffari 4399fe13df9SLeila Ghaffari // Create CEED Operator for each OutFlow faces 4409fe13df9SLeila Ghaffari if (nOut) { 4415603407fSLeila Ghaffari for (CeedInt i=0; i<nOut; i++) { 4425603407fSLeila Ghaffari ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, valueOut[i], numP_Sur, 4435603407fSLeila Ghaffari numQ_Sur, qdatasizeSur, &restrictqOut[i], &restrictxOut[i], 4445603407fSLeila Ghaffari &restrictqdiOut[i]); CHKERRQ(ierr); 4451e150236SLeila Ghaffari // Create the CEED vectors that will be needed in boundary setup 4465603407fSLeila Ghaffari CeedElemRestrictionGetNumElements(restrictqOut[i], &localNelemOut[i]); 4475603407fSLeila Ghaffari CeedVectorCreate(ceed, qdatasizeSur*localNelemOut[i]*NqptsSur, &qdataOut[i]); 4489fe13df9SLeila Ghaffari // Create the operator that builds the quadrature data for the OutFlow operator 4495603407fSLeila Ghaffari CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupOut[i]); 4505603407fSLeila Ghaffari CeedOperatorSetField(op_setupOut[i], "dx", restrictxOut[i], basisxSur, CEED_VECTOR_ACTIVE); 4515603407fSLeila Ghaffari CeedOperatorSetField(op_setupOut[i], "weight", CEED_ELEMRESTRICTION_NONE, 452ca3ac6ddSLeila Ghaffari basisxSur, CEED_VECTOR_NONE); 4535603407fSLeila Ghaffari CeedOperatorSetField(op_setupOut[i], "qdataSur", restrictqdiOut[i], 454ca3ac6ddSLeila Ghaffari CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 4559fe13df9SLeila Ghaffari // Create OutFlow operator 4565603407fSLeila Ghaffari CeedOperatorCreate(ceed, qf_applyOut, NULL, NULL, &op_applyOut[i]); 4575603407fSLeila Ghaffari CeedOperatorSetField(op_applyOut[i], "q", restrictqOut[i], basisqSur, CEED_VECTOR_ACTIVE); 4585603407fSLeila Ghaffari CeedOperatorSetField(op_applyOut[i], "qdataSur", restrictqdiOut[i], 4595603407fSLeila Ghaffari CEED_BASIS_COLLOCATED, qdataOut[i]); 4605603407fSLeila Ghaffari CeedOperatorSetField(op_applyOut[i], "x", restrictxOut[i], basisxSur, xcorners); 4615603407fSLeila Ghaffari CeedOperatorSetField(op_applyOut[i], "v", restrictqOut[i], basisqSur, CEED_VECTOR_ACTIVE); 4629fe13df9SLeila Ghaffari // Apply CEED operator for OutFlow setup 4635603407fSLeila Ghaffari CeedOperatorApply(op_setupOut[i], xcorners, qdataOut[i], CEED_REQUEST_IMMEDIATE); 4649fe13df9SLeila Ghaffari // Apply Sub-Operator for OutFlow BCs 4655603407fSLeila Ghaffari CeedCompositeOperatorAddSub(*op_apply, op_applyOut[i]); 466ca3ac6ddSLeila Ghaffari } 4679fe13df9SLeila Ghaffari } 4689fe13df9SLeila Ghaffari // Create CEED Operator for each InFlow faces 4699fe13df9SLeila Ghaffari if (nIn) { 4709fe13df9SLeila Ghaffari for (CeedInt i=0; i<nIn; i++) { 4719fe13df9SLeila Ghaffari ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, valueIn[i], numP_Sur, 4729fe13df9SLeila Ghaffari numQ_Sur, qdatasizeSur, &restrictqIn[i], &restrictxIn[i], 4739fe13df9SLeila Ghaffari &restrictqdiIn[i]); CHKERRQ(ierr); 4749fe13df9SLeila Ghaffari // Create the CEED vectors that will be needed in boundary setup 4759fe13df9SLeila Ghaffari CeedElemRestrictionGetNumElements(restrictqIn[i], &localNelemIn[i]); 4769fe13df9SLeila Ghaffari CeedVectorCreate(ceed, qdatasizeSur*localNelemIn[i]*NqptsSur, &qdataIn[i]); 4779fe13df9SLeila Ghaffari // Create the operator that builds the quadrature data for the InFlow operator 4789fe13df9SLeila Ghaffari CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupIn[i]); 4799fe13df9SLeila Ghaffari CeedOperatorSetField(op_setupIn[i], "dx", restrictxIn[i], basisxSur, CEED_VECTOR_ACTIVE); 4809fe13df9SLeila Ghaffari CeedOperatorSetField(op_setupIn[i], "weight", CEED_ELEMRESTRICTION_NONE, 4819fe13df9SLeila Ghaffari basisxSur, CEED_VECTOR_NONE); 4829fe13df9SLeila Ghaffari CeedOperatorSetField(op_setupIn[i], "qdataSur", restrictqdiIn[i], 4839fe13df9SLeila Ghaffari CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 4849fe13df9SLeila Ghaffari // Create InFlow operator 4859fe13df9SLeila Ghaffari CeedOperatorCreate(ceed, qf_applyIn, NULL, NULL, &op_applyIn[i]); 4869fe13df9SLeila Ghaffari CeedOperatorSetField(op_applyIn[i], "q", restrictqIn[i], basisqSur, CEED_VECTOR_ACTIVE); 4879fe13df9SLeila Ghaffari CeedOperatorSetField(op_applyIn[i], "qdataSur", restrictqdiIn[i], 4889fe13df9SLeila Ghaffari CEED_BASIS_COLLOCATED, qdataIn[i]); 4899fe13df9SLeila Ghaffari CeedOperatorSetField(op_applyIn[i], "x", restrictxIn[i], basisxSur, xcorners); 4909fe13df9SLeila Ghaffari CeedOperatorSetField(op_applyIn[i], "v", restrictqIn[i], basisqSur, CEED_VECTOR_ACTIVE); 4919fe13df9SLeila Ghaffari // Apply CEED operator for InFlow setup 4929fe13df9SLeila Ghaffari CeedOperatorApply(op_setupIn[i], xcorners, qdataIn[i], CEED_REQUEST_IMMEDIATE); 4939fe13df9SLeila Ghaffari // Apply Sub-Operator for InFlow BCs 4949fe13df9SLeila Ghaffari CeedCompositeOperatorAddSub(*op_apply, op_applyIn[i]); 4959fe13df9SLeila Ghaffari } 4969fe13df9SLeila Ghaffari } 4971e150236SLeila Ghaffari CeedVectorDestroy(&xcorners); 498ca3ac6ddSLeila Ghaffari } 499ca3ac6ddSLeila Ghaffari PetscFunctionReturn(0); 500ca3ac6ddSLeila Ghaffari } 501ca3ac6ddSLeila Ghaffari 502ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 503ccaff030SJeremy L Thompson PetscErrorCode ierr; 504ccaff030SJeremy L Thompson PetscInt m; 505ccaff030SJeremy L Thompson 506ccaff030SJeremy L Thompson PetscFunctionBeginUser; 507ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 508ccaff030SJeremy L Thompson ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 509ccaff030SJeremy L Thompson PetscFunctionReturn(0); 510ccaff030SJeremy L Thompson } 511ccaff030SJeremy L Thompson 512ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) { 513ccaff030SJeremy L Thompson PetscErrorCode ierr; 514ccaff030SJeremy L Thompson PetscInt mceed, mpetsc; 515ccaff030SJeremy L Thompson PetscScalar *a; 516ccaff030SJeremy L Thompson 517ccaff030SJeremy L Thompson PetscFunctionBeginUser; 518ccaff030SJeremy L Thompson ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 519ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 520ccaff030SJeremy L Thompson if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 52184d34d69SLeila Ghaffari "Cannot place PETSc Vec of length %D in CeedVector of length %D", 52284d34d69SLeila Ghaffari mpetsc, mceed); 523ccaff030SJeremy L Thompson ierr = VecGetArray(p, &a); CHKERRQ(ierr); 524ccaff030SJeremy L Thompson CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 525ccaff030SJeremy L Thompson PetscFunctionReturn(0); 526ccaff030SJeremy L Thompson } 527ccaff030SJeremy L Thompson 528ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 529ccaff030SJeremy L Thompson PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 530ccaff030SJeremy L Thompson Vec cellGeomFVM, Vec gradFVM) { 531ccaff030SJeremy L Thompson PetscErrorCode ierr; 532ccaff030SJeremy L Thompson Vec Qbc; 533ccaff030SJeremy L Thompson 534ccaff030SJeremy L Thompson PetscFunctionBegin; 535ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 536ccaff030SJeremy L Thompson ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 537ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 538ccaff030SJeremy L Thompson PetscFunctionReturn(0); 539ccaff030SJeremy L Thompson } 540ccaff030SJeremy L Thompson 541ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u) 542ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G 543ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 544ccaff030SJeremy L Thompson PetscErrorCode ierr; 545ccaff030SJeremy L Thompson User user = *(User *)userData; 546ccaff030SJeremy L Thompson PetscScalar *q, *g; 547ccaff030SJeremy L Thompson Vec Qloc, Gloc; 548ccaff030SJeremy L Thompson 549ccaff030SJeremy L Thompson // Global-to-local 550ccaff030SJeremy L Thompson PetscFunctionBeginUser; 551ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 552ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 553ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 554ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 555ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 556ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 557ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 558ccaff030SJeremy L Thompson 559ccaff030SJeremy L Thompson // Ceed Vectors 560ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 561ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 562ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 563ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 564ccaff030SJeremy L Thompson 565ccaff030SJeremy L Thompson // Apply CEED operator 566ccaff030SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 567ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 568ccaff030SJeremy L Thompson 569ccaff030SJeremy L Thompson // Restore vectors 570ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 571ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 572ccaff030SJeremy L Thompson 573ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 574ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 575ccaff030SJeremy L Thompson 576ccaff030SJeremy L Thompson // Inverse of the lumped mass matrix 577ccaff030SJeremy L Thompson ierr = VecPointwiseMult(G, G, user->M); // M is Minv 578ccaff030SJeremy L Thompson CHKERRQ(ierr); 579ccaff030SJeremy L Thompson 580ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 581ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 582ccaff030SJeremy L Thompson PetscFunctionReturn(0); 583ccaff030SJeremy L Thompson } 584ccaff030SJeremy L Thompson 585ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 586ccaff030SJeremy L Thompson void *userData) { 587ccaff030SJeremy L Thompson PetscErrorCode ierr; 588ccaff030SJeremy L Thompson User user = *(User *)userData; 589ccaff030SJeremy L Thompson const PetscScalar *q, *qdot; 590ccaff030SJeremy L Thompson PetscScalar *g; 591ccaff030SJeremy L Thompson Vec Qloc, Qdotloc, Gloc; 592ccaff030SJeremy L Thompson 593ccaff030SJeremy L Thompson // Global-to-local 594ccaff030SJeremy L Thompson PetscFunctionBeginUser; 595ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 596ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 597ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 598ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 599ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 600ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 601ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 602ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 603ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 604ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 605ccaff030SJeremy L Thompson 606ccaff030SJeremy L Thompson // Ceed Vectors 607ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 608ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 609ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 610ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 611ccaff030SJeremy L Thompson (PetscScalar *)q); 612ccaff030SJeremy L Thompson CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 613ccaff030SJeremy L Thompson (PetscScalar *)qdot); 614ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 615ccaff030SJeremy L Thompson 616ccaff030SJeremy L Thompson // Apply CEED operator 617ccaff030SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 618ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 619ccaff030SJeremy L Thompson 620ccaff030SJeremy L Thompson // Restore vectors 621ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 622ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 623ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 624ccaff030SJeremy L Thompson 625ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 626ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 627ccaff030SJeremy L Thompson 628ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 629ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 630ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 631ccaff030SJeremy L Thompson PetscFunctionReturn(0); 632ccaff030SJeremy L Thompson } 633ccaff030SJeremy L Thompson 634ccaff030SJeremy L Thompson // User provided TS Monitor 635ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 636ccaff030SJeremy L Thompson Vec Q, void *ctx) { 637ccaff030SJeremy L Thompson User user = ctx; 638ccaff030SJeremy L Thompson Vec Qloc; 639ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 640ccaff030SJeremy L Thompson PetscViewer viewer; 641ccaff030SJeremy L Thompson PetscErrorCode ierr; 642ccaff030SJeremy L Thompson 643ccaff030SJeremy L Thompson // Set up output 644ccaff030SJeremy L Thompson PetscFunctionBeginUser; 645ccaff030SJeremy L Thompson // Print every 'outputfreq' steps 646ccaff030SJeremy L Thompson if (stepno % user->outputfreq != 0) 647ccaff030SJeremy L Thompson PetscFunctionReturn(0); 648ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 649ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 650ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 651ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 652ccaff030SJeremy L Thompson 653ccaff030SJeremy L Thompson // Output 654ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 655ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 656ccaff030SJeremy L Thompson CHKERRQ(ierr); 657ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 658ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 659ccaff030SJeremy L Thompson ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 6609d801c56SJed Brown ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 661ccaff030SJeremy L Thompson if (user->dmviz) { 662ccaff030SJeremy L Thompson Vec Qrefined, Qrefined_loc; 663ccaff030SJeremy L Thompson char filepath_refined[PETSC_MAX_PATH_LEN]; 664ccaff030SJeremy L Thompson PetscViewer viewer_refined; 665ccaff030SJeremy L Thompson 666ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 667ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 668ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 669ccaff030SJeremy L Thompson CHKERRQ(ierr); 670ccaff030SJeremy L Thompson ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 671ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 672ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 673ccaff030SJeremy L Thompson CHKERRQ(ierr); 674ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 675ccaff030SJeremy L Thompson "%s/nsrefined-%03D.vtu", 676ccaff030SJeremy L Thompson user->outputfolder, stepno + user->contsteps); 677ccaff030SJeremy L Thompson CHKERRQ(ierr); 678ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 679ccaff030SJeremy L Thompson filepath_refined, 680ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 681ccaff030SJeremy L Thompson ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 682ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 683ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 684ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 685ccaff030SJeremy L Thompson } 686ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 687ccaff030SJeremy L Thompson 688ccaff030SJeremy L Thompson // Save data in a binary file for continuation of simulations 689ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 690ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 691ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 692ccaff030SJeremy L Thompson CHKERRQ(ierr); 693ccaff030SJeremy L Thompson ierr = VecView(Q, viewer); CHKERRQ(ierr); 694ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 695ccaff030SJeremy L Thompson 696ccaff030SJeremy L Thompson // Save time stamp 697ccaff030SJeremy L Thompson // Dimensionalize time back 698ccaff030SJeremy L Thompson time /= user->units->second; 699ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 700ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 701ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 702ccaff030SJeremy L Thompson CHKERRQ(ierr); 703ccaff030SJeremy L Thompson #if PETSC_VERSION_GE(3,13,0) 704ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 705ccaff030SJeremy L Thompson #else 706ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 707ccaff030SJeremy L Thompson #endif 708ccaff030SJeremy L Thompson CHKERRQ(ierr); 709ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 710ccaff030SJeremy L Thompson 711ccaff030SJeremy L Thompson PetscFunctionReturn(0); 712ccaff030SJeremy L Thompson } 713ccaff030SJeremy L Thompson 71484d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 715ccaff030SJeremy L Thompson CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 716ccaff030SJeremy L Thompson CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) { 717ccaff030SJeremy L Thompson PetscErrorCode ierr; 718ccaff030SJeremy L Thompson CeedVector multlvec; 719ccaff030SJeremy L Thompson Vec Multiplicity, MultiplicityLoc; 720ccaff030SJeremy L Thompson 721ccaff030SJeremy L Thompson ctxSetup->time = time; 722ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 723ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 724ccaff030SJeremy L Thompson CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 725ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 726ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 727ccaff030SJeremy L Thompson 728ccaff030SJeremy L Thompson // Fix multiplicity for output of ICs 729ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 730ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 731ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 732ccaff030SJeremy L Thompson CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 733ccaff030SJeremy L Thompson CeedVectorDestroy(&multlvec); 734ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 735ccaff030SJeremy L Thompson ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 736ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 737ccaff030SJeremy L Thompson CHKERRQ(ierr); 738ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 739ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 740ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 741ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 742ccaff030SJeremy L Thompson 743ccaff030SJeremy L Thompson PetscFunctionReturn(0); 744ccaff030SJeremy L Thompson } 745ccaff030SJeremy L Thompson 746ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 747ccaff030SJeremy L Thompson CeedElemRestriction restrictq, CeedBasis basisq, 748ccaff030SJeremy L Thompson CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 749ccaff030SJeremy L Thompson PetscErrorCode ierr; 750ccaff030SJeremy L Thompson CeedQFunction qf_mass; 751ccaff030SJeremy L Thompson CeedOperator op_mass; 752ccaff030SJeremy L Thompson CeedVector mceed; 753ccaff030SJeremy L Thompson Vec Mloc; 754ccaff030SJeremy L Thompson CeedInt ncompq, qdatasize; 755ccaff030SJeremy L Thompson 756ccaff030SJeremy L Thompson PetscFunctionBeginUser; 757ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 758ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 759ccaff030SJeremy L Thompson // Create the Q-function that defines the action of the mass operator 760ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 761ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 762ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 763ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 764ccaff030SJeremy L Thompson 765ccaff030SJeremy L Thompson // Create the mass operator 766ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 767ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 768ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", restrictqdi, 769ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 770ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 771ccaff030SJeremy L Thompson 772ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 773ccaff030SJeremy L Thompson ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 774ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 775ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 776ccaff030SJeremy L Thompson 777ccaff030SJeremy L Thompson { 778ccaff030SJeremy L Thompson // Compute a lumped mass matrix 779ccaff030SJeremy L Thompson CeedVector onesvec; 780ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 781ccaff030SJeremy L Thompson CeedVectorSetValue(onesvec, 1.0); 782ccaff030SJeremy L Thompson CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 783ccaff030SJeremy L Thompson CeedVectorDestroy(&onesvec); 784ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_mass); 785ccaff030SJeremy L Thompson CeedVectorDestroy(&mceed); 786ccaff030SJeremy L Thompson } 787ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_mass); 788ccaff030SJeremy L Thompson 789ccaff030SJeremy L Thompson ierr = VecZeroEntries(M); CHKERRQ(ierr); 790ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 791ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 792ccaff030SJeremy L Thompson 793ccaff030SJeremy L Thompson // Invert diagonally lumped mass vector for RHS function 794ccaff030SJeremy L Thompson ierr = VecReciprocal(M); CHKERRQ(ierr); 795ccaff030SJeremy L Thompson PetscFunctionReturn(0); 796ccaff030SJeremy L Thompson } 797ccaff030SJeremy L Thompson 79884d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 799ff6701fcSJed Brown SimpleBC bc, void *ctxSetup) { 800ccaff030SJeremy L Thompson PetscErrorCode ierr; 801ccaff030SJeremy L Thompson 802ccaff030SJeremy L Thompson PetscFunctionBeginUser; 803ccaff030SJeremy L Thompson { 804ccaff030SJeremy L Thompson // Configure the finite element space and boundary conditions 805ccaff030SJeremy L Thompson PetscFE fe; 806ccaff030SJeremy L Thompson PetscInt ncompq = 5; 807ff6701fcSJed Brown ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 808ff6701fcSJed Brown PETSC_FALSE, degree, PETSC_DECIDE, 80932ed2d11SJed Brown &fe); CHKERRQ(ierr); 810ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 811ccaff030SJeremy L Thompson ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 812ccaff030SJeremy L Thompson ierr = DMCreateDS(dm); CHKERRQ(ierr); 81307af6069Svaleriabarra { 81407af6069Svaleriabarra PetscInt comps[1] = {1}; 81507af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 81607af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[0], 81707af6069Svaleriabarra bc->slips[0], ctxSetup); CHKERRQ(ierr); 81807af6069Svaleriabarra comps[0] = 2; 81907af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 82007af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[1], 82107af6069Svaleriabarra bc->slips[1], ctxSetup); CHKERRQ(ierr); 82207af6069Svaleriabarra comps[0] = 3; 82307af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 82407af6069Svaleriabarra 1, comps, (void(*)(void))NULL, bc->nslip[2], 82507af6069Svaleriabarra bc->slips[2], ctxSetup); CHKERRQ(ierr); 82607af6069Svaleriabarra } 82784d34d69SLeila Ghaffari if (bc->userbc == PETSC_TRUE) { 82884d34d69SLeila Ghaffari for (PetscInt c = 0; c < 3; c++) { 82984d34d69SLeila Ghaffari for (PetscInt s = 0; s < bc->nslip[c]; s++) { 83084d34d69SLeila Ghaffari for (PetscInt w = 0; w < bc->nwall; w++) { 83184d34d69SLeila Ghaffari if (bc->slips[c][s] == bc->walls[w]) 83284d34d69SLeila Ghaffari SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 83384d34d69SLeila Ghaffari "Boundary condition already set on face %D!\n", bc->walls[w]); 83484d34d69SLeila Ghaffari 83584d34d69SLeila Ghaffari } 83684d34d69SLeila Ghaffari } 83784d34d69SLeila Ghaffari } 83884d34d69SLeila Ghaffari } 83984d34d69SLeila Ghaffari // Wall boundary conditions are zero energy density and zero flux for 84084d34d69SLeila Ghaffari // velocity in advection/advection2d, and zero velocity and zero flux 84184d34d69SLeila Ghaffari // for mass density and energy density in density_current 84284d34d69SLeila Ghaffari { 84384d34d69SLeila Ghaffari if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 84484d34d69SLeila Ghaffari PetscInt comps[1] = {4}; 84584d34d69SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 84684d34d69SLeila Ghaffari 1, comps, (void(*)(void))problem->bc, 84784d34d69SLeila Ghaffari bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 84884d34d69SLeila Ghaffari } else if (problem->bc == Exact_DC) { 84984d34d69SLeila Ghaffari PetscInt comps[3] = {1, 2, 3}; 85084d34d69SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 85184d34d69SLeila Ghaffari 3, comps, (void(*)(void))problem->bc, 85284d34d69SLeila Ghaffari bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 85384d34d69SLeila Ghaffari } else 85484d34d69SLeila Ghaffari SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, 85584d34d69SLeila Ghaffari "Undefined boundary conditions for this problem"); 85684d34d69SLeila Ghaffari } 857ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 858ccaff030SJeremy L Thompson CHKERRQ(ierr); 859ccaff030SJeremy L Thompson ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 860ccaff030SJeremy L Thompson } 861ccaff030SJeremy L Thompson { 862ccaff030SJeremy L Thompson // Empty name for conserved field (because there is only one field) 863ccaff030SJeremy L Thompson PetscSection section; 864ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 865ccaff030SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 866ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 867ccaff030SJeremy L Thompson CHKERRQ(ierr); 868ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 869ccaff030SJeremy L Thompson CHKERRQ(ierr); 870ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 871ccaff030SJeremy L Thompson CHKERRQ(ierr); 872ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 873ccaff030SJeremy L Thompson CHKERRQ(ierr); 874ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 875ccaff030SJeremy L Thompson CHKERRQ(ierr); 876ccaff030SJeremy L Thompson } 877ccaff030SJeremy L Thompson PetscFunctionReturn(0); 878ccaff030SJeremy L Thompson } 879ccaff030SJeremy L Thompson 880ccaff030SJeremy L Thompson int main(int argc, char **argv) { 881ccaff030SJeremy L Thompson PetscInt ierr; 882ccaff030SJeremy L Thompson MPI_Comm comm; 88384d34d69SLeila Ghaffari DM dm, dmcoord, dmviz; 884ccaff030SJeremy L Thompson Mat interpviz; 885ccaff030SJeremy L Thompson TS ts; 886ccaff030SJeremy L Thompson TSAdapt adapt; 887ccaff030SJeremy L Thompson User user; 888ccaff030SJeremy L Thompson Units units; 889ccaff030SJeremy L Thompson char ceedresource[4096] = "/cpu/self"; 89084d34d69SLeila Ghaffari PetscInt localNelemVol, lnodes, gnodes, steps; 891ccaff030SJeremy L Thompson const PetscInt ncompq = 5; 892ccaff030SJeremy L Thompson PetscMPIInt rank; 893ccaff030SJeremy L Thompson PetscScalar ftime; 894ccaff030SJeremy L Thompson Vec Q, Qloc, Xloc; 895ccaff030SJeremy L Thompson Ceed ceed; 89684d34d69SLeila Ghaffari CeedInt numP, numQ; 897cfa64770SLeila Ghaffari CeedVector xcorners, qdata, q0ceed; 89884d34d69SLeila Ghaffari CeedBasis basisx, basisxc, basisq; 89984d34d69SLeila Ghaffari CeedElemRestriction restrictx, restrictq, restrictqdi; 900cfa64770SLeila Ghaffari CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol; 901cfa64770SLeila Ghaffari CeedOperator op_setupVol, op_ics; 902ccaff030SJeremy L Thompson CeedScalar Rd; 90384d34d69SLeila Ghaffari CeedMemType memtyperequested; 904ccaff030SJeremy L Thompson PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 905ccaff030SJeremy L Thompson kgpersquaredms, Joulepercubicm; 906ccaff030SJeremy L Thompson problemType problemChoice; 907ccaff030SJeremy L Thompson problemData *problem = NULL; 9081184866aSLeila Ghaffari WindType wind_type; 909ccaff030SJeremy L Thompson StabilizationType stab; 91084d34d69SLeila Ghaffari testType testChoice; 91184d34d69SLeila Ghaffari testData *test = NULL; 91284d34d69SLeila Ghaffari PetscBool implicit; 913cb3e2689Svaleriabarra PetscInt viz_refine = 0; 914ccaff030SJeremy L Thompson struct SimpleBC_ bc = { 91584d34d69SLeila Ghaffari .nslip = {2, 2, 2}, 91684d34d69SLeila Ghaffari .slips = {{5, 6}, {3, 4}, {1, 2}} 917ccaff030SJeremy L Thompson }; 918ccaff030SJeremy L Thompson double start, cpu_time_used; 91984d34d69SLeila Ghaffari // Check PETSc CUDA support 92084d34d69SLeila Ghaffari PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 92184d34d69SLeila Ghaffari // *INDENT-OFF* 92284d34d69SLeila Ghaffari #ifdef PETSC_HAVE_CUDA 92384d34d69SLeila Ghaffari petschavecuda = PETSC_TRUE; 92484d34d69SLeila Ghaffari #else 92584d34d69SLeila Ghaffari petschavecuda = PETSC_FALSE; 92684d34d69SLeila Ghaffari #endif 92784d34d69SLeila Ghaffari // *INDENT-ON* 928ccaff030SJeremy L Thompson 929ccaff030SJeremy L Thompson // Create the libCEED contexts 930ccaff030SJeremy L Thompson PetscScalar meter = 1e-2; // 1 meter in scaled length units 931ccaff030SJeremy L Thompson PetscScalar second = 1e-2; // 1 second in scaled time units 932ccaff030SJeremy L Thompson PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 933ccaff030SJeremy L Thompson PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 934ccaff030SJeremy L Thompson CeedScalar theta0 = 300.; // K 935ccaff030SJeremy L Thompson CeedScalar thetaC = -15.; // K 936ccaff030SJeremy L Thompson CeedScalar P0 = 1.e5; // Pa 937*81f92cf0SLeila Ghaffari CeedScalar P_wind = 1.5e5; // Pa 938*81f92cf0SLeila Ghaffari CeedScalar rho_wind = 1.2; // Kg/m^3 939ccaff030SJeremy L Thompson CeedScalar N = 0.01; // 1/s 940ccaff030SJeremy L Thompson CeedScalar cv = 717.; // J/(kg K) 941ccaff030SJeremy L Thompson CeedScalar cp = 1004.; // J/(kg K) 942ccaff030SJeremy L Thompson CeedScalar g = 9.81; // m/s^2 943ccaff030SJeremy L Thompson CeedScalar lambda = -2./3.; // - 944ccaff030SJeremy L Thompson CeedScalar mu = 75.; // Pa s, dynamic viscosity 945ccaff030SJeremy L Thompson // mu = 75 is not physical for air, but is good for numerical stability 946ccaff030SJeremy L Thompson CeedScalar k = 0.02638; // W/(m K) 947ccaff030SJeremy L Thompson CeedScalar CtauS = 0.; // dimensionless 948ccaff030SJeremy L Thompson CeedScalar strong_form = 0.; // [0,1] 949ccaff030SJeremy L Thompson PetscScalar lx = 8000.; // m 950ccaff030SJeremy L Thompson PetscScalar ly = 8000.; // m 951ccaff030SJeremy L Thompson PetscScalar lz = 4000.; // m 952ccaff030SJeremy L Thompson CeedScalar rc = 1000.; // m (Radius of bubble) 953ccaff030SJeremy L Thompson PetscScalar resx = 1000.; // m (resolution in x) 954ccaff030SJeremy L Thompson PetscScalar resy = 1000.; // m (resolution in y) 955ccaff030SJeremy L Thompson PetscScalar resz = 1000.; // m (resolution in z) 956ccaff030SJeremy L Thompson PetscInt outputfreq = 10; // - 957ccaff030SJeremy L Thompson PetscInt contsteps = 0; // - 95884d34d69SLeila Ghaffari PetscInt degree = 1; // - 95984d34d69SLeila Ghaffari PetscInt qextra = 2; // - 960ea6e0f84SLeila Ghaffari PetscInt qextraSur = 2; // - 9611184866aSLeila Ghaffari PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0}; 962ccaff030SJeremy L Thompson 963ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 964ccaff030SJeremy L Thompson if (ierr) return ierr; 965ccaff030SJeremy L Thompson 966ccaff030SJeremy L Thompson // Allocate PETSc context 967ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 968ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 969ccaff030SJeremy L Thompson 970ccaff030SJeremy L Thompson // Parse command line options 971ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 972ccaff030SJeremy L Thompson ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 973ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 974ccaff030SJeremy L Thompson ierr = PetscOptionsString("-ceed", "CEED resource specifier", 975ccaff030SJeremy L Thompson NULL, ceedresource, ceedresource, 976ccaff030SJeremy L Thompson sizeof(ceedresource), NULL); CHKERRQ(ierr); 97784d34d69SLeila Ghaffari testChoice = TEST_NONE; 97884d34d69SLeila Ghaffari ierr = PetscOptionsEnum("-test", "Run tests", NULL, 97984d34d69SLeila Ghaffari testTypes, (PetscEnum)testChoice, 98084d34d69SLeila Ghaffari (PetscEnum *)&testChoice, 98184d34d69SLeila Ghaffari NULL); CHKERRQ(ierr); 98284d34d69SLeila Ghaffari test = &testOptions[testChoice]; 983ccaff030SJeremy L Thompson problemChoice = NS_DENSITY_CURRENT; 984ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 985ccaff030SJeremy L Thompson problemTypes, (PetscEnum)problemChoice, 986ccaff030SJeremy L Thompson (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 987ccaff030SJeremy L Thompson problem = &problemOptions[problemChoice]; 9881184866aSLeila Ghaffari ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection", 9891184866aSLeila Ghaffari NULL, WindTypes, (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION), 9901184866aSLeila Ghaffari (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr); 9911184866aSLeila Ghaffari PetscInt n = problem->dim; 99282c09b01SLeila Ghaffari PetscBool userWind; 9931184866aSLeila Ghaffari ierr = PetscOptionsRealArray("-problem_advection_wind_translation", "Constant wind vector", 99482c09b01SLeila Ghaffari NULL, wind, &n, &userWind); CHKERRQ(ierr); 99582c09b01SLeila Ghaffari if (wind_type == ADVECTION_WIND_ROTATION && userWind) { 99682c09b01SLeila Ghaffari ierr = PetscPrintf(comm, "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n"); 99782c09b01SLeila Ghaffari CHKERRQ(ierr); 9981184866aSLeila Ghaffari } 999ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 1000ccaff030SJeremy L Thompson StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 1001ccaff030SJeremy L Thompson (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 1002ccaff030SJeremy L Thompson ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 1003ccaff030SJeremy L Thompson NULL, implicit=PETSC_FALSE, &implicit, NULL); 1004ccaff030SJeremy L Thompson CHKERRQ(ierr); 100584d34d69SLeila Ghaffari if (!implicit && stab != STAB_NONE) { 100684d34d69SLeila Ghaffari ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 100784d34d69SLeila Ghaffari CHKERRQ(ierr); 100884d34d69SLeila Ghaffari } 1009ccaff030SJeremy L Thompson { 10107573aee6SJed Brown PetscInt len; 10117573aee6SJed Brown PetscBool flg; 1012ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray("-bc_wall", 1013ccaff030SJeremy L Thompson "Use wall boundary conditions on this list of faces", 1014ccaff030SJeremy L Thompson NULL, bc.walls, 1015ccaff030SJeremy L Thompson (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 1016ccaff030SJeremy L Thompson &len), &flg); CHKERRQ(ierr); 10177573aee6SJed Brown if (flg) { 10187573aee6SJed Brown bc.nwall = len; 10197573aee6SJed Brown // Using a no-slip wall disables automatic slip walls (they must be set explicitly) 10207573aee6SJed Brown bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0; 10217573aee6SJed Brown } 1022ccaff030SJeremy L Thompson for (PetscInt j=0; j<3; j++) { 1023ccaff030SJeremy L Thompson const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1024ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray(flags[j], 1025ccaff030SJeremy L Thompson "Use slip boundary conditions on this list of faces", 1026ccaff030SJeremy L Thompson NULL, bc.slips[j], 1027ccaff030SJeremy L Thompson (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 1028ccaff030SJeremy L Thompson &len), &flg); 1029ccaff030SJeremy L Thompson CHKERRQ(ierr); 103084d34d69SLeila Ghaffari if (flg) { 103184d34d69SLeila Ghaffari bc.nslip[j] = len; 103284d34d69SLeila Ghaffari bc.userbc = PETSC_TRUE; 103384d34d69SLeila Ghaffari } 1034ccaff030SJeremy L Thompson } 1035ccaff030SJeremy L Thompson } 1036cb3e2689Svaleriabarra ierr = PetscOptionsInt("-viz_refine", 1037cb3e2689Svaleriabarra "Regular refinement levels for visualization", 1038cb3e2689Svaleriabarra NULL, viz_refine, &viz_refine, NULL); 1039ccaff030SJeremy L Thompson CHKERRQ(ierr); 1040ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 1041ccaff030SJeremy L Thompson NULL, meter, &meter, NULL); CHKERRQ(ierr); 1042ccaff030SJeremy L Thompson meter = fabs(meter); 1043ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 1044ccaff030SJeremy L Thompson NULL, second, &second, NULL); CHKERRQ(ierr); 1045ccaff030SJeremy L Thompson second = fabs(second); 1046ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 1047ccaff030SJeremy L Thompson NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 1048ccaff030SJeremy L Thompson kilogram = fabs(kilogram); 1049ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_Kelvin", 1050ccaff030SJeremy L Thompson "1 Kelvin in scaled temperature units", 1051ccaff030SJeremy L Thompson NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 1052ccaff030SJeremy L Thompson Kelvin = fabs(Kelvin); 1053ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 1054ccaff030SJeremy L Thompson NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 1055ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 1056ccaff030SJeremy L Thompson NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 1057ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 1058ccaff030SJeremy L Thompson NULL, P0, &P0, NULL); CHKERRQ(ierr); 1059*81f92cf0SLeila Ghaffari ierr = PetscOptionsScalar("-P_wind", "Inflow wind pressure", 1060*81f92cf0SLeila Ghaffari NULL, P_wind, &P_wind, NULL); CHKERRQ(ierr); 1061*81f92cf0SLeila Ghaffari ierr = PetscOptionsScalar("-rho_wind", "Inflow wind density", 1062*81f92cf0SLeila Ghaffari NULL, rho_wind, &rho_wind, NULL); CHKERRQ(ierr); 1063ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 1064ccaff030SJeremy L Thompson NULL, N, &N, NULL); CHKERRQ(ierr); 1065ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 1066ccaff030SJeremy L Thompson NULL, cv, &cv, NULL); CHKERRQ(ierr); 1067ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 1068ccaff030SJeremy L Thompson NULL, cp, &cp, NULL); CHKERRQ(ierr); 1069ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 1070ccaff030SJeremy L Thompson NULL, g, &g, NULL); CHKERRQ(ierr); 1071ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lambda", 1072ccaff030SJeremy L Thompson "Stokes hypothesis second viscosity coefficient", 1073ccaff030SJeremy L Thompson NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 1074ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 1075ccaff030SJeremy L Thompson NULL, mu, &mu, NULL); CHKERRQ(ierr); 1076ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-k", "Thermal conductivity", 1077ccaff030SJeremy L Thompson NULL, k, &k, NULL); CHKERRQ(ierr); 1078ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-CtauS", 1079ccaff030SJeremy L Thompson "Scale coefficient for tau (nondimensional)", 1080ccaff030SJeremy L Thompson NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 108184d34d69SLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 108284d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 108384d34d69SLeila Ghaffari "Warning! Use -CtauS only with -stab su or -stab supg\n"); 108484d34d69SLeila Ghaffari CHKERRQ(ierr); 108584d34d69SLeila Ghaffari } 1086ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-strong_form", 1087ccaff030SJeremy L Thompson "Strong (1) or weak/integrated by parts (0) advection residual", 1088ccaff030SJeremy L Thompson NULL, strong_form, &strong_form, NULL); 1089ccaff030SJeremy L Thompson CHKERRQ(ierr); 109084d34d69SLeila Ghaffari if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) { 109184d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 109284d34d69SLeila Ghaffari "Warning! Problem density_current does not support -CtauS or -strong_form\n"); 109384d34d69SLeila Ghaffari CHKERRQ(ierr); 109484d34d69SLeila Ghaffari } 1095ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 1096ccaff030SJeremy L Thompson NULL, lx, &lx, NULL); CHKERRQ(ierr); 1097ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 1098ccaff030SJeremy L Thompson NULL, ly, &ly, NULL); CHKERRQ(ierr); 1099ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 1100ccaff030SJeremy L Thompson NULL, lz, &lz, NULL); CHKERRQ(ierr); 1101ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 1102ccaff030SJeremy L Thompson NULL, rc, &rc, NULL); CHKERRQ(ierr); 1103ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resx","Target resolution in x", 1104ccaff030SJeremy L Thompson NULL, resx, &resx, NULL); CHKERRQ(ierr); 1105ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resy","Target resolution in y", 1106ccaff030SJeremy L Thompson NULL, resy, &resy, NULL); CHKERRQ(ierr); 1107ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resz","Target resolution in z", 1108ccaff030SJeremy L Thompson NULL, resz, &resz, NULL); CHKERRQ(ierr); 110982c09b01SLeila Ghaffari n = problem->dim; 1110ccaff030SJeremy L Thompson center[0] = 0.5 * lx; 1111ccaff030SJeremy L Thompson center[1] = 0.5 * ly; 1112ccaff030SJeremy L Thompson center[2] = 0.5 * lz; 1113ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-center", "Location of bubble center", 1114ccaff030SJeremy L Thompson NULL, center, &n, NULL); CHKERRQ(ierr); 1115ccaff030SJeremy L Thompson n = problem->dim; 1116ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-dc_axis", 1117ccaff030SJeremy L Thompson "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 1118ccaff030SJeremy L Thompson NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 1119ccaff030SJeremy L Thompson { 1120ccaff030SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 1121ccaff030SJeremy L Thompson PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 1122ccaff030SJeremy L Thompson if (norm > 0) { 1123ccaff030SJeremy L Thompson for (int i=0; i<3; i++) dc_axis[i] /= norm; 1124ccaff030SJeremy L Thompson } 1125ccaff030SJeremy L Thompson } 1126ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-output_freq", 1127ccaff030SJeremy L Thompson "Frequency of output, in number of steps", 1128ccaff030SJeremy L Thompson NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 1129ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-continue", "Continue from previous solution", 1130ccaff030SJeremy L Thompson NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 113184d34d69SLeila Ghaffari ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 113284d34d69SLeila Ghaffari NULL, degree, °ree, NULL); CHKERRQ(ierr); 113384d34d69SLeila Ghaffari ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 113484d34d69SLeila Ghaffari NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 1135*81f92cf0SLeila Ghaffari PetscBool userQextraSur; 1136*81f92cf0SLeila Ghaffari ierr = PetscOptionsInt("-qextra_boundary", "Number of extra quadrature points on in/outflow faces", 1137*81f92cf0SLeila Ghaffari NULL, qextraSur, &qextraSur, &userQextraSur); CHKERRQ(ierr); 1138*81f92cf0SLeila Ghaffari if ((wind_type == ADVECTION_WIND_ROTATION || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) { 1139*81f92cf0SLeila Ghaffari ierr = PetscPrintf(comm, "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n"); 1140*81f92cf0SLeila Ghaffari CHKERRQ(ierr); 1141*81f92cf0SLeila Ghaffari } 114284d34d69SLeila Ghaffari ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr); 1143ccaff030SJeremy L Thompson ierr = PetscOptionsString("-of", "Output folder", 1144ccaff030SJeremy L Thompson NULL, user->outputfolder, user->outputfolder, 1145ccaff030SJeremy L Thompson sizeof(user->outputfolder), NULL); CHKERRQ(ierr); 114684d34d69SLeila Ghaffari memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 114784d34d69SLeila Ghaffari ierr = PetscOptionsEnum("-memtype", 114884d34d69SLeila Ghaffari "CEED MemType requested", NULL, 114984d34d69SLeila Ghaffari memTypes, (PetscEnum)memtyperequested, 115084d34d69SLeila Ghaffari (PetscEnum *)&memtyperequested, &setmemtyperequest); 115184d34d69SLeila Ghaffari CHKERRQ(ierr); 1152ccaff030SJeremy L Thompson ierr = PetscOptionsEnd(); CHKERRQ(ierr); 1153ccaff030SJeremy L Thompson 11541184866aSLeila Ghaffari // Setup BCs for Rotation or Translation wind types in Advection (3d) 11551184866aSLeila Ghaffari if (problemChoice == NS_ADVECTION) { 11561184866aSLeila Ghaffari switch (wind_type) { 11571184866aSLeila Ghaffari case ADVECTION_WIND_ROTATION: 11581184866aSLeila Ghaffari // No in/out-flow 11591184866aSLeila Ghaffari bc.ninflow = bc.noutflow = 0; 11601184866aSLeila Ghaffari break; 11611184866aSLeila Ghaffari case ADVECTION_WIND_TRANSLATION: 11621184866aSLeila Ghaffari // Face 6 is inflow and Face 5 is outflow 11631184866aSLeila Ghaffari bc.ninflow = bc.noutflow = 1; 11641184866aSLeila Ghaffari bc.inflow[0] = 6; bc.outflow[0] = 5; 11651184866aSLeila Ghaffari // Faces 3 and 4 are slip 11661184866aSLeila Ghaffari bc.nslip[0] = bc.nslip[2] = 0; bc.nslip[1] = 2; 11671184866aSLeila Ghaffari bc.slips[1][0] = 3; bc.slips[1][1] = 4; 11681184866aSLeila Ghaffari // Faces 1 and 2 are wall 11691184866aSLeila Ghaffari bc.nwall = 2; 11701184866aSLeila Ghaffari bc.walls[0] = 1; bc.walls[1] = 2; 11711184866aSLeila Ghaffari break; 11721184866aSLeila Ghaffari } 11731184866aSLeila Ghaffari } 11741184866aSLeila Ghaffari // Setup BCs for Rotation or Translation wind types in Advection (2d) 11751184866aSLeila Ghaffari if (problemChoice == NS_ADVECTION2D) { 11761184866aSLeila Ghaffari switch (wind_type) { 11771184866aSLeila Ghaffari case ADVECTION_WIND_ROTATION: 11781184866aSLeila Ghaffari break; 11791184866aSLeila Ghaffari case ADVECTION_WIND_TRANSLATION: 1180ed1c75c9SLeila Ghaffari bc.nwall = bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0; // No wall BCs, and slip BCs will be determined according to the wind vector. 1181ed1c75c9SLeila Ghaffari if (wind[0] == 0 && wind[1] == 0) { 1182ed1c75c9SLeila Ghaffari SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 1183ed1c75c9SLeila Ghaffari "No translation with problem_advection_wind_translation %f,%f. At least one of the elements needs to be non-zero.", 1184ed1c75c9SLeila Ghaffari wind[0], wind[1]); 1185ed1c75c9SLeila Ghaffari break;} else if (wind[0] == 0 && wind[1] > 0) { 1186ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 1; // Face 1 is inflow and Face 3 is outflow 1187ed1c75c9SLeila Ghaffari bc.inflow[0] = 1; bc.outflow[0] = 3; 1188ed1c75c9SLeila Ghaffari bc.nslip[0] = 2; // Faces 2 and 4 are slip 1189ed1c75c9SLeila Ghaffari bc.slips[0][0] = 2; bc.slips[0][1] = 4; 1190ed1c75c9SLeila Ghaffari break;} else if (wind[0] == 0 && wind[1] < 0) { 1191ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 1; // Face 3 is inflow and Face 1 is outflow 1192ed1c75c9SLeila Ghaffari bc.inflow[0] = 3; bc.outflow[0] = 1; 1193ed1c75c9SLeila Ghaffari bc.nslip[0] = 2; // Faces 2 and 4 are slip 1194ed1c75c9SLeila Ghaffari bc.slips[0][0] = 2; bc.slips[0][1] = 4; 1195ed1c75c9SLeila Ghaffari break;} else if (wind[0] > 0 && wind[1] == 0) { 1196ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 1; // Face 4 is inflow and Face 2 is outflow 11971184866aSLeila Ghaffari bc.inflow[0] = 4; bc.outflow[0] = 2; 1198ed1c75c9SLeila Ghaffari bc.nslip[1] = 2; // Faces 1 and 3 are slip 1199ed1c75c9SLeila Ghaffari bc.slips[1][0] = 1; bc.slips[1][1] = 3; 1200ed1c75c9SLeila Ghaffari break;} else if (wind[0] < 0 && wind[1] == 0) { 1201ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 1; // Face 2 is inflow and Face 4 is outflow 1202ed1c75c9SLeila Ghaffari bc.inflow[0] = 2; bc.outflow[0] = 4; 1203ed1c75c9SLeila Ghaffari bc.nslip[1] = 2; // Faces 1 and 3 are slip 1204ed1c75c9SLeila Ghaffari bc.slips[1][0] = 1; bc.slips[1][1] = 3; 1205ed1c75c9SLeila Ghaffari break;} else if (wind[0] > 0 && wind[1] > 0) { 1206ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 2; // Faces 1 and 4 are inflow and Faces 2 and 3 are outflow 1207ed1c75c9SLeila Ghaffari bc.inflow[0] = 1; bc.inflow[1] = 4; 1208ed1c75c9SLeila Ghaffari bc.outflow[0] = 2; bc.outflow[1] = 3; 1209ed1c75c9SLeila Ghaffari break;} else if (wind[0] > 0 && wind[1] < 0) { 1210ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 2; // Faces 3 and 4 are inflow and Faces 1 and 2 are outflow 1211ed1c75c9SLeila Ghaffari bc.inflow[0] = 3; bc.inflow[1] = 4; 1212ed1c75c9SLeila Ghaffari bc.outflow[0] = 1; bc.outflow[1] = 2; 1213ed1c75c9SLeila Ghaffari break;} else if (wind[0] < 0 && wind[1] > 0) { 1214ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 2; // Faces 1 and 2 are inflow and Faces 3 and 4 are outflow 1215ed1c75c9SLeila Ghaffari bc.inflow[0] = 1; bc.inflow[1] = 2; 1216ed1c75c9SLeila Ghaffari bc.outflow[0] = 3; bc.outflow[1] = 4; 1217ed1c75c9SLeila Ghaffari break;} else if (wind[0] < 0 && wind[1] < 0) { 1218ed1c75c9SLeila Ghaffari bc.ninflow = bc.noutflow = 2; // Faces 2 and 3 are inflow and Faces 1 and 4 are outflow 1219ed1c75c9SLeila Ghaffari bc.inflow[0] = 2; bc.inflow[1] = 3; 1220ed1c75c9SLeila Ghaffari bc.outflow[0] = 1; bc.outflow[1] = 4; 1221ed1c75c9SLeila Ghaffari break;} 12221184866aSLeila Ghaffari } 12231184866aSLeila Ghaffari } 12241184866aSLeila Ghaffari 1225ccaff030SJeremy L Thompson // Define derived units 1226ccaff030SJeremy L Thompson Pascal = kilogram / (meter * PetscSqr(second)); 1227ccaff030SJeremy L Thompson JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 1228ccaff030SJeremy L Thompson mpersquareds = meter / PetscSqr(second); 1229ccaff030SJeremy L Thompson WpermK = kilogram * meter / (pow(second,3) * Kelvin); 1230ccaff030SJeremy L Thompson kgpercubicm = kilogram / pow(meter,3); 1231ccaff030SJeremy L Thompson kgpersquaredms = kilogram / (PetscSqr(meter) * second); 1232ccaff030SJeremy L Thompson Joulepercubicm = kilogram / (meter * PetscSqr(second)); 1233ccaff030SJeremy L Thompson 1234ccaff030SJeremy L Thompson // Scale variables to desired units 1235ccaff030SJeremy L Thompson theta0 *= Kelvin; 1236ccaff030SJeremy L Thompson thetaC *= Kelvin; 1237ccaff030SJeremy L Thompson P0 *= Pascal; 1238*81f92cf0SLeila Ghaffari P_wind *= Pascal; 1239*81f92cf0SLeila Ghaffari rho_wind *= kgpercubicm; 1240ccaff030SJeremy L Thompson N *= (1./second); 1241ccaff030SJeremy L Thompson cv *= JperkgK; 1242ccaff030SJeremy L Thompson cp *= JperkgK; 1243ccaff030SJeremy L Thompson Rd = cp - cv; 1244ccaff030SJeremy L Thompson g *= mpersquareds; 1245ccaff030SJeremy L Thompson mu *= Pascal * second; 1246ccaff030SJeremy L Thompson k *= WpermK; 1247ccaff030SJeremy L Thompson lx = fabs(lx) * meter; 1248ccaff030SJeremy L Thompson ly = fabs(ly) * meter; 1249ccaff030SJeremy L Thompson lz = fabs(lz) * meter; 1250ccaff030SJeremy L Thompson rc = fabs(rc) * meter; 1251ccaff030SJeremy L Thompson resx = fabs(resx) * meter; 1252ccaff030SJeremy L Thompson resy = fabs(resy) * meter; 1253ccaff030SJeremy L Thompson resz = fabs(resz) * meter; 1254ccaff030SJeremy L Thompson for (int i=0; i<3; i++) center[i] *= meter; 1255ccaff030SJeremy L Thompson 1256ccaff030SJeremy L Thompson const CeedInt dim = problem->dim, ncompx = problem->dim, 1257cfa64770SLeila Ghaffari qdatasizeVol = problem->qdatasizeVol; 1258ccaff030SJeremy L Thompson // Set up the libCEED context 1259ccaff030SJeremy L Thompson struct SetupContext_ ctxSetup = { 1260ccaff030SJeremy L Thompson .theta0 = theta0, 1261ccaff030SJeremy L Thompson .thetaC = thetaC, 1262ccaff030SJeremy L Thompson .P0 = P0, 1263ccaff030SJeremy L Thompson .N = N, 1264ccaff030SJeremy L Thompson .cv = cv, 1265ccaff030SJeremy L Thompson .cp = cp, 1266ccaff030SJeremy L Thompson .Rd = Rd, 1267ccaff030SJeremy L Thompson .g = g, 1268ccaff030SJeremy L Thompson .rc = rc, 1269ccaff030SJeremy L Thompson .lx = lx, 1270ccaff030SJeremy L Thompson .ly = ly, 1271ccaff030SJeremy L Thompson .lz = lz, 1272ccaff030SJeremy L Thompson .center[0] = center[0], 1273ccaff030SJeremy L Thompson .center[1] = center[1], 1274ccaff030SJeremy L Thompson .center[2] = center[2], 1275ccaff030SJeremy L Thompson .dc_axis[0] = dc_axis[0], 1276ccaff030SJeremy L Thompson .dc_axis[1] = dc_axis[1], 1277ccaff030SJeremy L Thompson .dc_axis[2] = dc_axis[2], 12781184866aSLeila Ghaffari .wind[0] = wind[0], 12791184866aSLeila Ghaffari .wind[1] = wind[1], 12801184866aSLeila Ghaffari .wind[2] = wind[2], 1281ccaff030SJeremy L Thompson .time = 0, 12821184866aSLeila Ghaffari .wind_type = wind_type, 1283ccaff030SJeremy L Thompson }; 1284ccaff030SJeremy L Thompson 128584d34d69SLeila Ghaffari // Create the mesh 1286ccaff030SJeremy L Thompson { 1287ccaff030SJeremy L Thompson const PetscReal scale[3] = {lx, ly, lz}; 1288ccaff030SJeremy L Thompson ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 128984d34d69SLeila Ghaffari NULL, PETSC_TRUE, &dm); 1290ccaff030SJeremy L Thompson CHKERRQ(ierr); 1291ccaff030SJeremy L Thompson } 129284d34d69SLeila Ghaffari 129384d34d69SLeila Ghaffari // Distribute the mesh over processes 129484d34d69SLeila Ghaffari { 1295ccaff030SJeremy L Thompson DM dmDist = NULL; 1296ccaff030SJeremy L Thompson PetscPartitioner part; 1297ccaff030SJeremy L Thompson 1298ccaff030SJeremy L Thompson ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 1299ccaff030SJeremy L Thompson ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 1300ccaff030SJeremy L Thompson ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 1301ccaff030SJeremy L Thompson if (dmDist) { 1302ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1303ccaff030SJeremy L Thompson dm = dmDist; 1304ccaff030SJeremy L Thompson } 1305ccaff030SJeremy L Thompson } 1306ccaff030SJeremy L Thompson ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 1307ccaff030SJeremy L Thompson 130884d34d69SLeila Ghaffari // Setup DM 1309ccaff030SJeremy L Thompson ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 1310ccaff030SJeremy L Thompson ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 131184d34d69SLeila Ghaffari ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr); 131284d34d69SLeila Ghaffari 131384d34d69SLeila Ghaffari // Refine DM for high-order viz 1314ccaff030SJeremy L Thompson dmviz = NULL; 1315ccaff030SJeremy L Thompson interpviz = NULL; 1316ccaff030SJeremy L Thompson if (viz_refine) { 1317ff6701fcSJed Brown DM dmhierarchy[viz_refine+1]; 1318ff6701fcSJed Brown 1319ccaff030SJeremy L Thompson ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 1320ff6701fcSJed Brown dmhierarchy[0] = dm; 132184d34d69SLeila Ghaffari for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 1322ff6701fcSJed Brown Mat interp_next; 1323ff6701fcSJed Brown 1324ff6701fcSJed Brown ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 1325ccaff030SJeremy L Thompson CHKERRQ(ierr); 1326ff6701fcSJed Brown ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 1327ff6701fcSJed Brown d = (d + 1) / 2; 1328ff6701fcSJed Brown if (i + 1 == viz_refine) d = 1; 1329ff6701fcSJed Brown ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr); 1330ff6701fcSJed Brown ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 1331ff6701fcSJed Brown &interp_next, NULL); CHKERRQ(ierr); 1332ff6701fcSJed Brown if (!i) interpviz = interp_next; 1333ff6701fcSJed Brown else { 1334ff6701fcSJed Brown Mat C; 1335ff6701fcSJed Brown ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 1336ff6701fcSJed Brown PETSC_DECIDE, &C); CHKERRQ(ierr); 1337ff6701fcSJed Brown ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 1338ff6701fcSJed Brown ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1339ff6701fcSJed Brown interpviz = C; 1340ff6701fcSJed Brown } 1341ff6701fcSJed Brown } 1342cb3e2689Svaleriabarra for (PetscInt i=1; i<viz_refine; i++) { 1343ff6701fcSJed Brown ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 1344cb3e2689Svaleriabarra } 1345ff6701fcSJed Brown dmviz = dmhierarchy[viz_refine]; 1346ccaff030SJeremy L Thompson } 1347ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 1348ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 1349ccaff030SJeremy L Thompson ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 1350ccaff030SJeremy L Thompson lnodes /= ncompq; 1351ccaff030SJeremy L Thompson 135284d34d69SLeila Ghaffari // Initialize CEED 135384d34d69SLeila Ghaffari CeedInit(ceedresource, &ceed); 135484d34d69SLeila Ghaffari // Set memtype 135584d34d69SLeila Ghaffari CeedMemType memtypebackend; 135684d34d69SLeila Ghaffari CeedGetPreferredMemType(ceed, &memtypebackend); 135784d34d69SLeila Ghaffari // Check memtype compatibility 135884d34d69SLeila Ghaffari if (!setmemtyperequest) 135984d34d69SLeila Ghaffari memtyperequested = memtypebackend; 136084d34d69SLeila Ghaffari else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 136184d34d69SLeila Ghaffari SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 136284d34d69SLeila Ghaffari "PETSc was not built with CUDA. " 136384d34d69SLeila Ghaffari "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 136484d34d69SLeila Ghaffari 136584d34d69SLeila Ghaffari // Set number of 1D nodes and quadrature points 136684d34d69SLeila Ghaffari numP = degree + 1; 136784d34d69SLeila Ghaffari numQ = numP + qextra; 136884d34d69SLeila Ghaffari 136984d34d69SLeila Ghaffari // Print summary 137084d34d69SLeila Ghaffari if (testChoice == TEST_NONE) { 1371ccaff030SJeremy L Thompson CeedInt gdofs, odofs; 1372ccaff030SJeremy L Thompson int comm_size; 1373ccaff030SJeremy L Thompson char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 1374ccaff030SJeremy L Thompson ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 1375ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 137684d34d69SLeila Ghaffari gnodes = gdofs/ncompq; 1377ccaff030SJeremy L Thompson ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 1378ccaff030SJeremy L Thompson ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 1379ccaff030SJeremy L Thompson sizeof(box_faces_str), NULL); CHKERRQ(ierr); 138084d34d69SLeila Ghaffari const char *usedresource; 138184d34d69SLeila Ghaffari CeedGetResource(ceed, &usedresource); 1382ccaff030SJeremy L Thompson 138384d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 138484d34d69SLeila Ghaffari "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 138584d34d69SLeila Ghaffari " rank(s) : %d\n" 138684d34d69SLeila Ghaffari " Problem:\n" 138784d34d69SLeila Ghaffari " Problem Name : %s\n" 138884d34d69SLeila Ghaffari " Stabilization : %s\n" 138984d34d69SLeila Ghaffari " PETSc:\n" 139084d34d69SLeila Ghaffari " Box Faces : %s\n" 139184d34d69SLeila Ghaffari " libCEED:\n" 139284d34d69SLeila Ghaffari " libCEED Backend : %s\n" 139384d34d69SLeila Ghaffari " libCEED Backend MemType : %s\n" 139484d34d69SLeila Ghaffari " libCEED User Requested MemType : %s\n" 139584d34d69SLeila Ghaffari " Mesh:\n" 139684d34d69SLeila Ghaffari " Number of 1D Basis Nodes (P) : %d\n" 139784d34d69SLeila Ghaffari " Number of 1D Quadrature Points (Q) : %d\n" 139884d34d69SLeila Ghaffari " Global DoFs : %D\n" 139984d34d69SLeila Ghaffari " Owned DoFs : %D\n" 140084d34d69SLeila Ghaffari " DoFs per node : %D\n" 140184d34d69SLeila Ghaffari " Global nodes : %D\n" 140284d34d69SLeila Ghaffari " Owned nodes : %D\n", 140384d34d69SLeila Ghaffari comm_size, problemTypes[problemChoice], 140484d34d69SLeila Ghaffari StabilizationTypes[stab], box_faces_str, usedresource, 140584d34d69SLeila Ghaffari CeedMemTypes[memtypebackend], 140684d34d69SLeila Ghaffari (setmemtyperequest) ? 140784d34d69SLeila Ghaffari CeedMemTypes[memtyperequested] : "none", 140884d34d69SLeila Ghaffari numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes); 140984d34d69SLeila Ghaffari CHKERRQ(ierr); 14100c6c0b13SLeila Ghaffari } 14110c6c0b13SLeila Ghaffari 1412ccaff030SJeremy L Thompson // Set up global mass vector 1413ccaff030SJeremy L Thompson ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 1414ccaff030SJeremy L Thompson 141584d34d69SLeila Ghaffari // Set up libCEED 1416ccaff030SJeremy L Thompson // CEED Bases 1417ccaff030SJeremy L Thompson CeedInit(ceedresource, &ceed); 141884d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 141984d34d69SLeila Ghaffari &basisq); 142084d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 142184d34d69SLeila Ghaffari &basisx); 142284d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 142384d34d69SLeila Ghaffari CEED_GAUSS_LOBATTO, &basisxc); 1424ccaff030SJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 1425ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 1426ccaff030SJeremy L Thompson CHKERRQ(ierr); 1427ccaff030SJeremy L Thompson 1428ccaff030SJeremy L Thompson // CEED Restrictions 14291e150236SLeila Ghaffari ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ, 143084d34d69SLeila Ghaffari qdatasizeVol, &restrictq, &restrictx, 143184d34d69SLeila Ghaffari &restrictqdi); CHKERRQ(ierr); 1432ccaff030SJeremy L Thompson 1433ccaff030SJeremy L Thompson ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1434ccaff030SJeremy L Thompson ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1435ccaff030SJeremy L Thompson 1436ccaff030SJeremy L Thompson // Create the CEED vectors that will be needed in setup 1437bd910870SLeila Ghaffari CeedInt NqptsVol; 143884d34d69SLeila Ghaffari CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol); 143984d34d69SLeila Ghaffari CeedElemRestrictionGetNumElements(restrictq, &localNelemVol); 14408b982baeSLeila Ghaffari CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata); 144184d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1442ccaff030SJeremy L Thompson 1443ccaff030SJeremy L Thompson // Create the Q-function that builds the quadrature data for the NS operator 1444ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc, 1445ea6e0f84SLeila Ghaffari &qf_setupVol); 1446ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD); 1447ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT); 14488b982baeSLeila Ghaffari CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1449ccaff030SJeremy L Thompson 1450ccaff030SJeremy L Thompson // Create the Q-function that sets the ICs of the operator 1451ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1452ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1453ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1454ccaff030SJeremy L Thompson 1455ea6e0f84SLeila Ghaffari qf_rhsVol = NULL; 1456ea6e0f84SLeila Ghaffari if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator 1457ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs, 1458ea6e0f84SLeila Ghaffari problem->applyVol_rhs_loc, &qf_rhsVol); 1459ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP); 1460ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 14618b982baeSLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1462ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP); 1463ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP); 1464ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1465ccaff030SJeremy L Thompson } 1466ccaff030SJeremy L Thompson 1467ea6e0f84SLeila Ghaffari qf_ifunctionVol = NULL; 1468ea6e0f84SLeila Ghaffari if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction 1469ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction, 1470ea6e0f84SLeila Ghaffari problem->applyVol_ifunction_loc, &qf_ifunctionVol); 1471ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP); 1472ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 1473ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP); 14748b982baeSLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1475ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP); 1476ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP); 1477ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1478ccaff030SJeremy L Thompson } 1479ccaff030SJeremy L Thompson 1480ccaff030SJeremy L Thompson // Create the operator that builds the quadrature data for the NS operator 1481ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol); 148284d34d69SLeila Ghaffari CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1483ea6e0f84SLeila Ghaffari CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE, 148484d34d69SLeila Ghaffari basisx, CEED_VECTOR_NONE); 148584d34d69SLeila Ghaffari CeedOperatorSetField(op_setupVol, "qdata", restrictqdi, 1486ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1487ccaff030SJeremy L Thompson 1488ccaff030SJeremy L Thompson // Create the operator that sets the ICs 1489ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 149084d34d69SLeila Ghaffari CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 149184d34d69SLeila Ghaffari CeedOperatorSetField(op_ics, "q0", restrictq, 1492ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1493ccaff030SJeremy L Thompson 149484d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 149584d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 149684d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1497ccaff030SJeremy L Thompson 1498ea6e0f84SLeila Ghaffari if (qf_rhsVol) { // Create the RHS physics operator 1499ccaff030SJeremy L Thompson CeedOperator op; 1500ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op); 150184d34d69SLeila Ghaffari CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 150284d34d69SLeila Ghaffari CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 150384d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdata", restrictqdi, 15048b982baeSLeila Ghaffari CEED_BASIS_COLLOCATED, qdata); 150584d34d69SLeila Ghaffari CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 150684d34d69SLeila Ghaffari CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 150784d34d69SLeila Ghaffari CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1508d3630711SJed Brown user->op_rhs_vol = op; 1509ccaff030SJeremy L Thompson } 1510ccaff030SJeremy L Thompson 1511ea6e0f84SLeila Ghaffari if (qf_ifunctionVol) { // Create the IFunction operator 1512ccaff030SJeremy L Thompson CeedOperator op; 1513ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op); 151484d34d69SLeila Ghaffari CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 151584d34d69SLeila Ghaffari CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 151684d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 151784d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdata", restrictqdi, 15188b982baeSLeila Ghaffari CEED_BASIS_COLLOCATED, qdata); 151984d34d69SLeila Ghaffari CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 152084d34d69SLeila Ghaffari CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 152184d34d69SLeila Ghaffari CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1522d3630711SJed Brown user->op_ifunction_vol = op; 1523ccaff030SJeremy L Thompson } 1524ccaff030SJeremy L Thompson 1525cfa64770SLeila Ghaffari //--------------------------------------------------------------------------------------// 15269fe13df9SLeila Ghaffari // In/OutFlow Boundary Conditions 15276a0edaf9SLeila Ghaffari //--------------------------------------------------------------------------------------// 15286a0edaf9SLeila Ghaffari // Set up CEED for the boundaries 15296a0edaf9SLeila Ghaffari CeedInt height = 1; 15306a0edaf9SLeila Ghaffari CeedInt dimSur = dim - height; 15311e150236SLeila Ghaffari CeedInt numP_Sur = degree + 1; 15321e150236SLeila Ghaffari CeedInt numQ_Sur = numP_Sur + qextraSur; 1533cfa64770SLeila Ghaffari const CeedInt qdatasizeSur = problem->qdatasizeSur; 1534cfa64770SLeila Ghaffari CeedBasis basisxSur, basisxcSur, basisqSur; 1535cfa64770SLeila Ghaffari CeedInt NqptsSur; 15369fe13df9SLeila Ghaffari CeedQFunction qf_setupSur, qf_rhsOut, qf_ifunctionOut, qf_rhsIn, qf_ifunctionIn; 1537cfa64770SLeila Ghaffari 1538cfa64770SLeila Ghaffari // CEED bases for the boundaries 15396a0edaf9SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS, 15406a0edaf9SLeila Ghaffari &basisqSur); 15416a0edaf9SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS, 15426a0edaf9SLeila Ghaffari &basisxSur); 15436a0edaf9SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur, 15446a0edaf9SLeila Ghaffari CEED_GAUSS_LOBATTO, &basisxcSur); 15456a0edaf9SLeila Ghaffari CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur); 15466a0edaf9SLeila Ghaffari 1547cfa64770SLeila Ghaffari // Create the Q-function that builds the quadrature data for the Surface operator 15486a0edaf9SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc, 15496a0edaf9SLeila Ghaffari &qf_setupSur); 15506a0edaf9SLeila Ghaffari CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD); 15516a0edaf9SLeila Ghaffari CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT); 15526a0edaf9SLeila Ghaffari CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 15536a0edaf9SLeila Ghaffari 15549fe13df9SLeila Ghaffari // Creat Q-Function for OutFlow BCs 15555603407fSLeila Ghaffari qf_rhsOut = NULL; 15569fe13df9SLeila Ghaffari if (problem->applyOut_rhs) { // Create the Q-function that defines the action of the RHS operator 15575603407fSLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyOut_rhs, 15585603407fSLeila Ghaffari problem->applyOut_rhs_loc, &qf_rhsOut); 15595603407fSLeila Ghaffari CeedQFunctionAddInput(qf_rhsOut, "q", ncompq, CEED_EVAL_INTERP); 15605603407fSLeila Ghaffari CeedQFunctionAddInput(qf_rhsOut, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 15615603407fSLeila Ghaffari CeedQFunctionAddInput(qf_rhsOut, "x", ncompx, CEED_EVAL_INTERP); 15625603407fSLeila Ghaffari CeedQFunctionAddOutput(qf_rhsOut, "v", ncompq, CEED_EVAL_INTERP); 15636a0edaf9SLeila Ghaffari } 15645603407fSLeila Ghaffari qf_ifunctionOut = NULL; 15655603407fSLeila Ghaffari if (problem->applyOut_ifunction) { // Create the Q-function that defines the action of the IFunction 15665603407fSLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyOut_ifunction, 15675603407fSLeila Ghaffari problem->applyOut_ifunction_loc, &qf_ifunctionOut); 15685603407fSLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionOut, "q", ncompq, CEED_EVAL_INTERP); 15695603407fSLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionOut, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 15705603407fSLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionOut, "x", ncompx, CEED_EVAL_INTERP); 15715603407fSLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionOut, "v", ncompq, CEED_EVAL_INTERP); 15726a0edaf9SLeila Ghaffari } 15739fe13df9SLeila Ghaffari 15749fe13df9SLeila Ghaffari // Creat Q-Function for InFlow BCs 15759fe13df9SLeila Ghaffari qf_rhsIn = NULL; 15769fe13df9SLeila Ghaffari if (problem->applyIn_rhs) { // Create the Q-function that defines the action of the RHS operator 15779fe13df9SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyIn_rhs, 15789fe13df9SLeila Ghaffari problem->applyIn_rhs_loc, &qf_rhsIn); 15799fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_rhsIn, "q", ncompq, CEED_EVAL_INTERP); 15809fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_rhsIn, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 15819fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_rhsIn, "x", ncompx, CEED_EVAL_INTERP); 15829fe13df9SLeila Ghaffari CeedQFunctionAddOutput(qf_rhsIn, "v", ncompq, CEED_EVAL_INTERP); 15839fe13df9SLeila Ghaffari } 15849fe13df9SLeila Ghaffari qf_ifunctionIn = NULL; 15859fe13df9SLeila Ghaffari if (problem->applyIn_ifunction) { // Create the Q-function that defines the action of the IFunction 15869fe13df9SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyIn_ifunction, 15879fe13df9SLeila Ghaffari problem->applyIn_ifunction_loc, &qf_ifunctionIn); 15889fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionIn, "q", ncompq, CEED_EVAL_INTERP); 15899fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionIn, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 15909fe13df9SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionIn, "x", ncompx, CEED_EVAL_INTERP); 15919fe13df9SLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionIn, "v", ncompq, CEED_EVAL_INTERP); 15929fe13df9SLeila Ghaffari } 15939fe13df9SLeila Ghaffari 15949fe13df9SLeila Ghaffari // Create CEED Operator for the whole domain 15959fe13df9SLeila Ghaffari if (!implicit) 15969fe13df9SLeila Ghaffari ierr = CreateOperatorForDomain(ceed, dm, user->op_rhs_vol, qf_rhsOut, qf_rhsIn, qf_setupSur, height, 15979fe13df9SLeila Ghaffari bc.noutflow, bc.outflow, bc.ninflow, bc.inflow, numP_Sur, numQ_Sur, qdatasizeSur, 15981e150236SLeila Ghaffari NqptsSur, basisxSur, basisqSur, &user->op_rhs); 1599ca3ac6ddSLeila Ghaffari CHKERRQ(ierr); 16009fe13df9SLeila Ghaffari if (implicit) 16019fe13df9SLeila Ghaffari ierr = CreateOperatorForDomain(ceed, dm, user->op_ifunction_vol, qf_ifunctionOut, qf_ifunctionIn, qf_setupSur, height, 16029fe13df9SLeila Ghaffari bc.noutflow, bc.outflow, bc.ninflow, bc.inflow, numP_Sur, numQ_Sur, qdatasizeSur, 16031e150236SLeila Ghaffari NqptsSur, basisxSur, basisqSur, &user->op_ifunction); 1604ca3ac6ddSLeila Ghaffari CHKERRQ(ierr); 16059fe13df9SLeila Ghaffari 1606cfa64770SLeila Ghaffari //--------------------------------------------------------------------------------------// 1607ccaff030SJeremy L Thompson CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup); 1608ccaff030SJeremy L Thompson CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd}; 1609*81f92cf0SLeila Ghaffari CeedScalar ctxIn[5] = {cv, cp, Rd, P_wind, rho_wind}; 1610ccaff030SJeremy L Thompson struct Advection2dContext_ ctxAdvection2d = { 1611ccaff030SJeremy L Thompson .CtauS = CtauS, 1612ccaff030SJeremy L Thompson .strong_form = strong_form, 1613ccaff030SJeremy L Thompson .stabilization = stab, 1614ccaff030SJeremy L Thompson }; 1615ccaff030SJeremy L Thompson switch (problemChoice) { 1616ccaff030SJeremy L Thompson case NS_DENSITY_CURRENT: 1617ea6e0f84SLeila Ghaffari if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS); 1618ea6e0f84SLeila Ghaffari if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS, 1619ccaff030SJeremy L Thompson sizeof ctxNS); 1620ccaff030SJeremy L Thompson break; 1621ccaff030SJeremy L Thompson case NS_ADVECTION: 1622ccaff030SJeremy L Thompson case NS_ADVECTION2D: 1623ea6e0f84SLeila Ghaffari if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d, 1624ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 1625ea6e0f84SLeila Ghaffari if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d, 1626ccaff030SJeremy L Thompson sizeof ctxAdvection2d); 16275603407fSLeila Ghaffari if (qf_rhsOut) CeedQFunctionSetContext(qf_rhsOut, &ctxAdvection2d, 16286a0edaf9SLeila Ghaffari sizeof ctxAdvection2d); 16295603407fSLeila Ghaffari if (qf_ifunctionOut) CeedQFunctionSetContext(qf_ifunctionOut, &ctxAdvection2d, 16306a0edaf9SLeila Ghaffari sizeof ctxAdvection2d); 16319fe13df9SLeila Ghaffari if (qf_rhsIn) CeedQFunctionSetContext(qf_rhsIn, &ctxIn, sizeof ctxIn); 16329fe13df9SLeila Ghaffari if (qf_ifunctionIn) CeedQFunctionSetContext(qf_ifunctionIn, &ctxIn, sizeof ctxIn); 1633ccaff030SJeremy L Thompson } 1634ccaff030SJeremy L Thompson 1635ccaff030SJeremy L Thompson // Set up PETSc context 1636ccaff030SJeremy L Thompson // Set up units structure 1637ccaff030SJeremy L Thompson units->meter = meter; 1638ccaff030SJeremy L Thompson units->kilogram = kilogram; 1639ccaff030SJeremy L Thompson units->second = second; 1640ccaff030SJeremy L Thompson units->Kelvin = Kelvin; 1641ccaff030SJeremy L Thompson units->Pascal = Pascal; 1642ccaff030SJeremy L Thompson units->JperkgK = JperkgK; 1643ccaff030SJeremy L Thompson units->mpersquareds = mpersquareds; 1644ccaff030SJeremy L Thompson units->WpermK = WpermK; 1645ccaff030SJeremy L Thompson units->kgpercubicm = kgpercubicm; 1646ccaff030SJeremy L Thompson units->kgpersquaredms = kgpersquaredms; 1647ccaff030SJeremy L Thompson units->Joulepercubicm = Joulepercubicm; 1648ccaff030SJeremy L Thompson 1649ccaff030SJeremy L Thompson // Set up user structure 1650ccaff030SJeremy L Thompson user->comm = comm; 1651ccaff030SJeremy L Thompson user->outputfreq = outputfreq; 1652ccaff030SJeremy L Thompson user->contsteps = contsteps; 1653ccaff030SJeremy L Thompson user->units = units; 1654ccaff030SJeremy L Thompson user->dm = dm; 1655ccaff030SJeremy L Thompson user->dmviz = dmviz; 1656ccaff030SJeremy L Thompson user->interpviz = interpviz; 1657ccaff030SJeremy L Thompson user->ceed = ceed; 1658ccaff030SJeremy L Thompson 16598b982baeSLeila Ghaffari // Calculate qdata and ICs 1660ccaff030SJeremy L Thompson // Set up state global and local vectors 1661ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1662ccaff030SJeremy L Thompson 1663cfa64770SLeila Ghaffari ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1664ccaff030SJeremy L Thompson 1665ccaff030SJeremy L Thompson // Apply Setup Ceed Operators 1666ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 16678b982baeSLeila Ghaffari CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 166884d34d69SLeila Ghaffari ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1669ccaff030SJeremy L Thompson user->M); CHKERRQ(ierr); 1670ccaff030SJeremy L Thompson 167184d34d69SLeila Ghaffari ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 167284d34d69SLeila Ghaffari &ctxSetup, 0.0); CHKERRQ(ierr); 1673ccaff030SJeremy L Thompson if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1674ccaff030SJeremy L Thompson // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1675ccaff030SJeremy L Thompson // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1676ccaff030SJeremy L Thompson // slower execution. 1677ccaff030SJeremy L Thompson Vec Qbc; 1678ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1679ccaff030SJeremy L Thompson ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1680ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1681ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1682ccaff030SJeremy L Thompson ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1683ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1684ccaff030SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)dm, 168584d34d69SLeila Ghaffari "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 168684d34d69SLeila Ghaffari CHKERRQ(ierr); 1687ccaff030SJeremy L Thompson } 1688ccaff030SJeremy L Thompson 1689ccaff030SJeremy L Thompson MPI_Comm_rank(comm, &rank); 1690ccaff030SJeremy L Thompson if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);} 1691ccaff030SJeremy L Thompson // Gather initial Q values 1692ccaff030SJeremy L Thompson // In case of continuation of simulation, set up initial values from binary file 1693ccaff030SJeremy L Thompson if (contsteps) { // continue from existent solution 1694ccaff030SJeremy L Thompson PetscViewer viewer; 1695ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1696ccaff030SJeremy L Thompson // Read input 1697ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1698ccaff030SJeremy L Thompson user->outputfolder); 1699ccaff030SJeremy L Thompson CHKERRQ(ierr); 1700ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1701ccaff030SJeremy L Thompson CHKERRQ(ierr); 1702ccaff030SJeremy L Thompson ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1703ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1704ccaff030SJeremy L Thompson } 1705ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1706ccaff030SJeremy L Thompson 1707ccaff030SJeremy L Thompson // Create and setup TS 1708ccaff030SJeremy L Thompson ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1709ccaff030SJeremy L Thompson ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1710ccaff030SJeremy L Thompson if (implicit) { 1711ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1712ccaff030SJeremy L Thompson if (user->op_ifunction) { 1713ccaff030SJeremy L Thompson ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1714ccaff030SJeremy L Thompson } else { // Implicit integrators can fall back to using an RHSFunction 1715ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1716ccaff030SJeremy L Thompson } 1717ccaff030SJeremy L Thompson } else { 1718ccaff030SJeremy L Thompson if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, 1719ccaff030SJeremy L Thompson "Problem does not provide RHSFunction"); 1720ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1721ccaff030SJeremy L Thompson ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1722ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1723ccaff030SJeremy L Thompson } 1724ccaff030SJeremy L Thompson ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1725ccaff030SJeremy L Thompson ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1726ccaff030SJeremy L Thompson ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 172784d34d69SLeila Ghaffari if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);} 1728ccaff030SJeremy L Thompson ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1729ccaff030SJeremy L Thompson ierr = TSAdaptSetStepLimits(adapt, 1730ccaff030SJeremy L Thompson 1.e-12 * units->second, 1731ccaff030SJeremy L Thompson 1.e2 * units->second); CHKERRQ(ierr); 1732ccaff030SJeremy L Thompson ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1733ccaff030SJeremy L Thompson if (!contsteps) { // print initial condition 173484d34d69SLeila Ghaffari if (testChoice == TEST_NONE) { 1735ccaff030SJeremy L Thompson ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1736ccaff030SJeremy L Thompson } 1737ccaff030SJeremy L Thompson } else { // continue from time of last output 1738ccaff030SJeremy L Thompson PetscReal time; 1739ccaff030SJeremy L Thompson PetscInt count; 1740ccaff030SJeremy L Thompson PetscViewer viewer; 1741ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1742ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1743ccaff030SJeremy L Thompson user->outputfolder); CHKERRQ(ierr); 1744ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1745ccaff030SJeremy L Thompson CHKERRQ(ierr); 1746ccaff030SJeremy L Thompson ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1747ccaff030SJeremy L Thompson CHKERRQ(ierr); 1748ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1749ccaff030SJeremy L Thompson ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1750ccaff030SJeremy L Thompson } 175184d34d69SLeila Ghaffari if (testChoice == TEST_NONE) { 1752ccaff030SJeremy L Thompson ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1753ccaff030SJeremy L Thompson } 1754ccaff030SJeremy L Thompson 1755ccaff030SJeremy L Thompson // Solve 1756ccaff030SJeremy L Thompson start = MPI_Wtime(); 1757ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1758ccaff030SJeremy L Thompson ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1759ccaff030SJeremy L Thompson cpu_time_used = MPI_Wtime() - start; 1760ccaff030SJeremy L Thompson ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1761ccaff030SJeremy L Thompson ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1762ccaff030SJeremy L Thompson comm); CHKERRQ(ierr); 176384d34d69SLeila Ghaffari if (testChoice == TEST_NONE) { 1764ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 176584d34d69SLeila Ghaffari "Time taken for solution (sec): %g\n", 1766ccaff030SJeremy L Thompson (double)cpu_time_used); CHKERRQ(ierr); 1767ccaff030SJeremy L Thompson } 1768ccaff030SJeremy L Thompson 1769ccaff030SJeremy L Thompson // Get error 177084d34d69SLeila Ghaffari if (problem->non_zero_time && testChoice == TEST_NONE) { 1771ccaff030SJeremy L Thompson Vec Qexact, Qexactloc; 1772ccaff030SJeremy L Thompson PetscReal norm; 1773ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1774ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1775ccaff030SJeremy L Thompson ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1776ccaff030SJeremy L Thompson 177784d34d69SLeila Ghaffari ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 177884d34d69SLeila Ghaffari restrictq, &ctxSetup, ftime); CHKERRQ(ierr); 1779ccaff030SJeremy L Thompson 1780ccaff030SJeremy L Thompson ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1781ccaff030SJeremy L Thompson ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1782cfa64770SLeila Ghaffari CeedVectorDestroy(&q0ceed); 1783ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1784ccaff030SJeremy L Thompson "Max Error: %g\n", 1785ccaff030SJeremy L Thompson (double)norm); CHKERRQ(ierr); 178684d34d69SLeila Ghaffari // Clean up vectors 178784d34d69SLeila Ghaffari ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 178884d34d69SLeila Ghaffari ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1789ccaff030SJeremy L Thompson } 1790ccaff030SJeremy L Thompson 1791ccaff030SJeremy L Thompson // Output Statistics 1792ccaff030SJeremy L Thompson ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr); 179384d34d69SLeila Ghaffari if (testChoice == TEST_NONE) { 1794ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1795ccaff030SJeremy L Thompson "Time integrator took %D time steps to reach final time %g\n", 1796ccaff030SJeremy L Thompson steps, (double)ftime); CHKERRQ(ierr); 1797ccaff030SJeremy L Thompson } 179884d34d69SLeila Ghaffari // Output numerical values from command line 179984d34d69SLeila Ghaffari ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr); 180084d34d69SLeila Ghaffari 180184d34d69SLeila Ghaffari // compare reference solution values with current run 180284d34d69SLeila Ghaffari if (testChoice != TEST_NONE) { 180384d34d69SLeila Ghaffari PetscViewer viewer; 180484d34d69SLeila Ghaffari // Read reference file 180584d34d69SLeila Ghaffari Vec Qref; 180684d34d69SLeila Ghaffari PetscReal error, Qrefnorm; 180784d34d69SLeila Ghaffari ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr); 180884d34d69SLeila Ghaffari ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer); 180984d34d69SLeila Ghaffari CHKERRQ(ierr); 181084d34d69SLeila Ghaffari ierr = VecLoad(Qref, viewer); CHKERRQ(ierr); 181184d34d69SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 181284d34d69SLeila Ghaffari 181384d34d69SLeila Ghaffari // Compute error with respect to reference solution 181484d34d69SLeila Ghaffari ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr); 181584d34d69SLeila Ghaffari ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr); 181684d34d69SLeila Ghaffari ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr); 181784d34d69SLeila Ghaffari ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr); 181884d34d69SLeila Ghaffari ierr = VecDestroy(&Qref); CHKERRQ(ierr); 181984d34d69SLeila Ghaffari // Check error 182084d34d69SLeila Ghaffari if (error > test->testtol) { 182184d34d69SLeila Ghaffari ierr = PetscPrintf(PETSC_COMM_WORLD, 182284d34d69SLeila Ghaffari "Test failed with error norm %g\n", 182384d34d69SLeila Ghaffari (double)error); CHKERRQ(ierr); 182484d34d69SLeila Ghaffari } 182584d34d69SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 182684d34d69SLeila Ghaffari } 18279cf88b28Svaleriabarra 1828ccaff030SJeremy L Thompson // Clean up libCEED 18298b982baeSLeila Ghaffari CeedVectorDestroy(&qdata); 1830ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qceed); 1831ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qdotceed); 1832ccaff030SJeremy L Thompson CeedVectorDestroy(&user->gceed); 1833ccaff030SJeremy L Thompson CeedVectorDestroy(&xcorners); 183484d34d69SLeila Ghaffari CeedBasisDestroy(&basisq); 183584d34d69SLeila Ghaffari CeedBasisDestroy(&basisx); 183684d34d69SLeila Ghaffari CeedBasisDestroy(&basisxc); 183784d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictq); 183884d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictx); 183984d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictqdi); 1840ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_setupVol); 1841ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ics); 1842ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_rhsVol); 1843ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_ifunctionVol); 1844ea6e0f84SLeila Ghaffari CeedOperatorDestroy(&op_setupVol); 1845ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_ics); 18466a0edaf9SLeila Ghaffari CeedOperatorDestroy(&user->op_rhs_vol); 18476a0edaf9SLeila Ghaffari CeedOperatorDestroy(&user->op_ifunction_vol); 18486a0edaf9SLeila Ghaffari CeedDestroy(&ceed); 18496a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisqSur); 18506a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisxSur); 18516a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisxcSur); 18526a0edaf9SLeila Ghaffari CeedQFunctionDestroy(&qf_setupSur); 18535603407fSLeila Ghaffari CeedQFunctionDestroy(&qf_rhsOut); 18545603407fSLeila Ghaffari CeedQFunctionDestroy(&qf_ifunctionOut); 18559fe13df9SLeila Ghaffari CeedQFunctionDestroy(&qf_rhsIn); 18569fe13df9SLeila Ghaffari CeedQFunctionDestroy(&qf_ifunctionIn); 1857ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_rhs); 1858ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_ifunction); 1859ccaff030SJeremy L Thompson 1860ccaff030SJeremy L Thompson // Clean up PETSc 1861ccaff030SJeremy L Thompson ierr = VecDestroy(&Q); CHKERRQ(ierr); 1862ccaff030SJeremy L Thompson ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1863ccaff030SJeremy L Thompson ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1864ccaff030SJeremy L Thompson ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1865ccaff030SJeremy L Thompson ierr = TSDestroy(&ts); CHKERRQ(ierr); 1866ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1867ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 1868ccaff030SJeremy L Thompson ierr = PetscFree(user); CHKERRQ(ierr); 1869ccaff030SJeremy L Thompson return PetscFinalize(); 1870ccaff030SJeremy L Thompson } 1871ccaff030SJeremy L Thompson 1872