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