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