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/cuda -problem advection -degree 1 33 // 34 //TESTARGS(name="dc_explicit") -ceed {ceed_resource} -test -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-explicit.bin 35 //TESTARGS(name="dc_implicit_stab_none") -ceed {ceed_resource} -test -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-implicit-stab-none.bin 36 //TESTARGS(name="adv_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection -strong_form 1 -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-explicit-strong.bin 37 //TESTARGS(name="adv_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -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 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-implicit-stab-supg.bin 38 //TESTARGS(name="adv_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -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 -problem_advection_wind translation -problem_advection_wind_translation .53,-1.33,-2.65 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-translation-implicit-stab-su.bin 39 //TESTARGS(name="adv2d_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection2d -strong_form 1 -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-explicit-strong.bin 40 //TESTARGS(name="adv2d_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -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 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-implicit-stab-supg.bin 41 //TESTARGS(name="adv2d_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -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 -problem_advection_wind translation -problem_advection_wind_translation .53,-1.33,0 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-translation-implicit-stab-su.bin 42 43 /// @file 44 /// Navier-Stokes example using PETSc 45 46 const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 47 48 #include <ceed.h> 49 #include <petscdmplex.h> 50 #include <petscsys.h> 51 #include <petscts.h> 52 #include <stdbool.h> 53 #include "advection.h" 54 #include "advection2d.h" 55 #include "common.h" 56 #include "euler-vortex.h" 57 #include "densitycurrent.h" 58 #include "setup-boundary.h" 59 60 #if PETSC_VERSION_LT(3,14,0) 61 # define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i) 62 # define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i) 63 #endif 64 65 #if PETSC_VERSION_LT(3,14,0) 66 # define DMAddBoundary(a,b,c,d,e,f,g,h,i,j,k,l) DMAddBoundary(a,b,c,d,e,f,g,h,j,k,l) 67 #endif 68 69 // MemType Options 70 static const char *const memTypes[] = { 71 "host", 72 "device", 73 "memType", "CEED_MEM_", NULL 74 }; 75 76 // Problem Options 77 typedef enum { 78 NS_DENSITY_CURRENT = 0, 79 NS_ADVECTION = 1, 80 NS_ADVECTION2D = 2, 81 NS_EULER_VORTEX = 3, 82 } problemType; 83 static const char *const problemTypes[] = { 84 "density_current", 85 "advection", 86 "advection2d", 87 "euler_vortex", 88 "problemType", "NS_", NULL 89 }; 90 91 // Wind Options for Advection 92 typedef enum { 93 ADVECTION_WIND_ROTATION = 0, 94 ADVECTION_WIND_TRANSLATION = 1, 95 } WindType; 96 static const char *const WindTypes[] = { 97 "rotation", 98 "translation", 99 "WindType", "ADVECTION_WIND_", NULL 100 }; 101 102 typedef enum { 103 STAB_NONE = 0, 104 STAB_SU = 1, // Streamline Upwind 105 STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 106 } StabilizationType; 107 static const char *const StabilizationTypes[] = { 108 "none", 109 "SU", 110 "SUPG", 111 "StabilizationType", "STAB_", NULL 112 }; 113 114 // Problem specific data 115 typedef struct { 116 CeedInt dim, qdatasizeVol, qdatasizeSur; 117 CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction, 118 applySur; 119 PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 120 PetscScalar[], void *); 121 const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc, 122 *applyVol_ifunction_loc, *applySur_loc; 123 const bool non_zero_time; 124 } problemData; 125 126 problemData problemOptions[] = { 127 [NS_DENSITY_CURRENT] = { 128 .dim = 3, 129 .qdatasizeVol = 10, 130 .qdatasizeSur = 4, 131 .setupVol = Setup, 132 .setupVol_loc = Setup_loc, 133 .setupSur = SetupBoundary, 134 .setupSur_loc = SetupBoundary_loc, 135 .ics = ICsDC, 136 .ics_loc = ICsDC_loc, 137 .applyVol_rhs = DC, 138 .applyVol_rhs_loc = DC_loc, 139 .applyVol_ifunction = IFunction_DC, 140 .applyVol_ifunction_loc = IFunction_DC_loc, 141 .bc = Exact_DC, 142 .non_zero_time = PETSC_FALSE, 143 }, 144 [NS_ADVECTION] = { 145 .dim = 3, 146 .qdatasizeVol = 10, 147 .qdatasizeSur = 4, 148 .setupVol = Setup, 149 .setupVol_loc = Setup_loc, 150 .setupSur = SetupBoundary, 151 .setupSur_loc = SetupBoundary_loc, 152 .ics = ICsAdvection, 153 .ics_loc = ICsAdvection_loc, 154 .applyVol_rhs = Advection, 155 .applyVol_rhs_loc = Advection_loc, 156 .applyVol_ifunction = IFunction_Advection, 157 .applyVol_ifunction_loc = IFunction_Advection_loc, 158 .applySur = Advection_Sur, 159 .applySur_loc = Advection_Sur_loc, 160 .bc = Exact_Advection, 161 .non_zero_time = PETSC_FALSE, 162 }, 163 [NS_ADVECTION2D] = { 164 .dim = 2, 165 .qdatasizeVol = 5, 166 .qdatasizeSur = 3, 167 .setupVol = Setup2d, 168 .setupVol_loc = Setup2d_loc, 169 .setupSur = SetupBoundary2d, 170 .setupSur_loc = SetupBoundary2d_loc, 171 .ics = ICsAdvection2d, 172 .ics_loc = ICsAdvection2d_loc, 173 .applyVol_rhs = Advection2d, 174 .applyVol_rhs_loc = Advection2d_loc, 175 .applyVol_ifunction = IFunction_Advection2d, 176 .applyVol_ifunction_loc = IFunction_Advection2d_loc, 177 .applySur = Advection2d_Sur, 178 .applySur_loc = Advection2d_Sur_loc, 179 .bc = Exact_Advection2d, 180 .non_zero_time = PETSC_TRUE, 181 }, 182 [NS_EULER_VORTEX] = { 183 .dim = 3, 184 .qdatasizeVol = 10, 185 .qdatasizeSur = 4, 186 .setupVol = Setup, 187 .setupVol_loc = Setup_loc, 188 .setupSur = SetupBoundary, 189 .setupSur_loc = SetupBoundary_loc, 190 .ics = ICsEuler, 191 .ics_loc = ICsEuler_loc, 192 .applyVol_rhs = Euler, 193 .applyVol_rhs_loc = Euler_loc, 194 .bc = Exact_Euler, 195 .non_zero_time = PETSC_TRUE, 196 }, 197 }; 198 199 // PETSc user data 200 typedef struct User_ *User; 201 typedef struct Units_ *Units; 202 203 struct User_ { 204 MPI_Comm comm; 205 PetscInt outputfreq; 206 DM dm; 207 DM dmviz; 208 Mat interpviz; 209 Ceed ceed; 210 Units units; 211 CeedVector qceed, qdotceed, gceed; 212 CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction; 213 Vec M; 214 char outputdir[PETSC_MAX_PATH_LEN]; 215 PetscInt contsteps; 216 EulerContext ctxEulerData; 217 }; 218 219 struct Units_ { 220 // fundamental units 221 PetscScalar meter; 222 PetscScalar kilogram; 223 PetscScalar second; 224 PetscScalar Kelvin; 225 // derived units 226 PetscScalar Pascal; 227 PetscScalar JperkgK; 228 PetscScalar mpersquareds; 229 PetscScalar WpermK; 230 PetscScalar kgpercubicm; 231 PetscScalar kgpersquaredms; 232 PetscScalar Joulepercubicm; 233 PetscScalar Joule; 234 }; 235 236 typedef struct SimpleBC_ *SimpleBC; 237 struct SimpleBC_ { 238 PetscInt nwall, nslip[3]; 239 PetscInt walls[6], slips[3][6]; 240 PetscBool userbc; 241 }; 242 243 // Essential BC dofs are encoded in closure indices as -(i+1). 244 static PetscInt Involute(PetscInt i) { 245 return i >= 0 ? i : -(i+1); 246 } 247 248 // Utility function to create local CEED restriction 249 static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 250 CeedInt height, DMLabel domainLabel, CeedInt value, 251 CeedElemRestriction *Erestrict) { 252 253 PetscSection section; 254 PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth; 255 DMLabel depthLabel; 256 IS depthIS, iterIS; 257 Vec Uloc; 258 const PetscInt *iterIndices; 259 PetscErrorCode ierr; 260 261 PetscFunctionBeginUser; 262 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 263 dim -= height; 264 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 265 ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 266 PetscInt ncomp[nfields], fieldoff[nfields+1]; 267 fieldoff[0] = 0; 268 for (PetscInt f=0; f<nfields; f++) { 269 ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 270 fieldoff[f+1] = fieldoff[f] + ncomp[f]; 271 } 272 273 ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 274 ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr); 275 ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr); 276 if (domainLabel) { 277 IS domainIS; 278 ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr); 279 if (domainIS) { // domainIS is non-empty 280 ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr); 281 ierr = ISDestroy(&domainIS); CHKERRQ(ierr); 282 } else { // domainIS is NULL (empty) 283 iterIS = NULL; 284 } 285 ierr = ISDestroy(&depthIS); CHKERRQ(ierr); 286 } else { 287 iterIS = depthIS; 288 } 289 if (iterIS) { 290 ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr); 291 ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr); 292 } else { 293 Nelem = 0; 294 iterIndices = NULL; 295 } 296 ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 297 for (p=0,eoffset=0; p<Nelem; p++) { 298 PetscInt c = iterIndices[p]; 299 PetscInt numindices, *indices, nnodes; 300 ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 301 &numindices, &indices, NULL, NULL); 302 CHKERRQ(ierr); 303 bool flip = false; 304 if (height > 0) { 305 PetscInt numCells, numFaces, start = -1; 306 const PetscInt *orients, *faces, *cells; 307 ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 308 ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr); 309 if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 310 "Expected one cell in support of exterior face, but got %D cells", 311 numCells); 312 ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 313 ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr); 314 for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;} 315 if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, 316 "Could not find face %D in cone of its support", 317 c); 318 ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 319 if (orients[start] < 0) flip = true; 320 } 321 if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 322 PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 323 c); 324 nnodes = numindices / fieldoff[nfields]; 325 for (PetscInt i=0; i<nnodes; i++) { 326 PetscInt ii = i; 327 if (flip) { 328 if (P == nnodes) ii = nnodes - 1 - i; 329 else if (P*P == nnodes) { 330 PetscInt row = i / P, col = i % P; 331 ii = row + col * P; 332 } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, 333 "No support for flipping point with %D nodes != P (%D) or P^2", 334 nnodes, P); 335 } 336 // Check that indices are blocked by node and thus can be coalesced as a single field with 337 // fieldoff[nfields] = sum(ncomp) components. 338 for (PetscInt f=0; f<nfields; f++) { 339 for (PetscInt j=0; j<ncomp[f]; j++) { 340 if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j]) 341 != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j) 342 SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 343 "Cell %D closure indices not interlaced for node %D field %D component %D", 344 c, ii, f, j); 345 } 346 } 347 // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 348 PetscInt loc = Involute(indices[ii*ncomp[0]]); 349 erestrict[eoffset++] = loc; 350 } 351 ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 352 &numindices, &indices, NULL, NULL); 353 CHKERRQ(ierr); 354 } 355 if (eoffset != Nelem*PetscPowInt(P, dim)) 356 SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 357 "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 358 PetscPowInt(P, dim),eoffset); 359 if (iterIS) { 360 ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr); 361 } 362 ierr = ISDestroy(&iterIS); CHKERRQ(ierr); 363 364 ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 365 ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 366 ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 367 CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields], 368 1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 369 Erestrict); 370 ierr = PetscFree(erestrict); CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 374 // Utility function to get Ceed Restriction for each domain 375 static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, 376 DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize, 377 CeedElemRestriction *restrictq, CeedElemRestriction *restrictx, 378 CeedElemRestriction *restrictqdi) { 379 380 DM dmcoord; 381 CeedInt dim, localNelem; 382 CeedInt Qdim; 383 PetscErrorCode ierr; 384 385 PetscFunctionBeginUser; 386 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 387 dim -= height; 388 Qdim = CeedIntPow(Q, dim); 389 ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 390 ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 391 CHKERRQ(ierr); 392 ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, 393 restrictq); 394 CHKERRQ(ierr); 395 ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, 396 restrictx); 397 CHKERRQ(ierr); 398 CeedElemRestrictionGetNumElements(*restrictq, &localNelem); 399 CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim, 400 qdatasize, qdatasize*localNelem*Qdim, 401 CEED_STRIDES_BACKEND, restrictqdi); 402 PetscFunctionReturn(0); 403 } 404 405 // Utility function to create CEED Composite Operator for the entire domain 406 static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, 407 problemType problemChoice, WindType wind_type, CeedOperator op_applyVol, 408 CeedQFunction qf_applySur, CeedQFunction qf_setupSur,CeedInt height, 409 CeedInt numP_Sur, CeedInt numQ_Sur, CeedInt qdatasizeSur, CeedInt NqptsSur, 410 CeedBasis basisxSur, CeedBasis basisqSur, CeedOperator *op_apply) { 411 412 CeedInt dim, nFace; 413 PetscInt lsize; 414 Vec Xloc; 415 CeedVector xcorners; 416 DMLabel domainLabel; 417 PetscScalar *x; 418 PetscErrorCode ierr; 419 420 PetscFunctionBeginUser; 421 // Composite Operaters 422 CeedCompositeOperatorCreate(ceed, op_apply); 423 // --Apply a Sub-Operator for the volume 424 CeedCompositeOperatorAddSub(*op_apply, op_applyVol); 425 426 // Required data for in/outflow BCs 427 ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 428 ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr); 429 ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr); 430 ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr); 431 CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x); 432 ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr); 433 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 434 435 if (wind_type == ADVECTION_WIND_TRANSLATION) { 436 // Ignore wall and slip BCs 437 bc->nwall = 0; 438 bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0; 439 440 // Set number of faces 441 if (dim == 2) nFace = 4; 442 if (dim == 3) nFace = 6; 443 444 // Create CEED Operator for each boundary face 445 PetscInt localNelemSur[6]; 446 CeedVector qdataSur[6]; 447 CeedOperator op_setupSur[6], op_applySur[6]; 448 CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6]; 449 450 for (CeedInt i=0; i<nFace; i++) { 451 ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur, 452 numQ_Sur, qdatasizeSur, &restrictqSur[i], 453 &restrictxSur[i], &restrictqdiSur[i]); 454 CHKERRQ(ierr); 455 // Create the CEED vectors that will be needed in Boundary setup 456 CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]); 457 CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, 458 &qdataSur[i]); 459 // Create the operator that builds the quadrature data for the Boundary operator 460 CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]); 461 CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, 462 CEED_VECTOR_ACTIVE); 463 CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE, 464 basisxSur, CEED_VECTOR_NONE); 465 CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i], 466 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 467 // Create Boundary operator 468 CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]); 469 CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur, 470 CEED_VECTOR_ACTIVE); 471 CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i], 472 CEED_BASIS_COLLOCATED, qdataSur[i]); 473 CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, 474 xcorners); 475 CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur, 476 CEED_VECTOR_ACTIVE); 477 // Apply CEED operator for Boundary setup 478 CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i], 479 CEED_REQUEST_IMMEDIATE); 480 // --Apply Sub-Operator for the Boundary 481 CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]); 482 } 483 CeedVectorDestroy(&xcorners); 484 } 485 PetscFunctionReturn(0); 486 } 487 488 static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 489 PetscErrorCode ierr; 490 PetscInt m; 491 492 PetscFunctionBeginUser; 493 ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 494 ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497 498 static int VectorPlacePetscVec(CeedVector c, Vec p) { 499 PetscErrorCode ierr; 500 PetscInt mceed, mpetsc; 501 PetscScalar *a; 502 503 PetscFunctionBeginUser; 504 ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 505 ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 506 if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 507 "Cannot place PETSc Vec of length %D in CeedVector of length %D", 508 mpetsc, mceed); 509 ierr = VecGetArray(p, &a); CHKERRQ(ierr); 510 CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 511 PetscFunctionReturn(0); 512 } 513 514 static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 515 PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 516 Vec cellGeomFVM, Vec gradFVM) { 517 PetscErrorCode ierr; 518 Vec Qbc; 519 520 PetscFunctionBegin; 521 ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 522 ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 523 ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 524 PetscFunctionReturn(0); 525 } 526 527 // This is the RHS of the ODE, given as u_t = G(t,u) 528 // This function takes in a state vector Q and writes into G 529 static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 530 PetscErrorCode ierr; 531 User user = *(User *)userData; 532 PetscScalar *q, *g; 533 Vec Qloc, Gloc; 534 535 // Global-to-local 536 PetscFunctionBeginUser; 537 user->ctxEulerData->currentTime = t; 538 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 539 ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 540 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 541 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 542 ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 543 NULL, NULL, NULL); CHKERRQ(ierr); 544 ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 545 546 // Ceed Vectors 547 ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 548 ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 549 CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 550 CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 551 552 // Apply CEED operator 553 CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 554 CEED_REQUEST_IMMEDIATE); 555 556 // Restore vectors 557 ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 558 ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 559 560 ierr = VecZeroEntries(G); CHKERRQ(ierr); 561 ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 562 563 // Inverse of the lumped mass matrix 564 ierr = VecPointwiseMult(G, G, user->M); // M is Minv 565 CHKERRQ(ierr); 566 567 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 568 ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 569 PetscFunctionReturn(0); 570 } 571 572 static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 573 void *userData) { 574 PetscErrorCode ierr; 575 User user = *(User *)userData; 576 const PetscScalar *q, *qdot; 577 PetscScalar *g; 578 Vec Qloc, Qdotloc, Gloc; 579 580 // Global-to-local 581 PetscFunctionBeginUser; 582 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 583 ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 584 ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 585 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 586 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 587 ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 588 NULL, NULL, NULL); CHKERRQ(ierr); 589 ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 590 ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 591 ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 592 593 // Ceed Vectors 594 ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 595 ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 596 ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 597 CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 598 (PetscScalar *)q); 599 CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 600 (PetscScalar *)qdot); 601 CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 602 603 // Apply CEED operator 604 CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 605 CEED_REQUEST_IMMEDIATE); 606 607 // Restore vectors 608 ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 609 ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 610 ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 611 612 ierr = VecZeroEntries(G); CHKERRQ(ierr); 613 ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 614 615 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 616 ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 617 ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 // User provided TS Monitor 622 static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 623 Vec Q, void *ctx) { 624 User user = ctx; 625 Vec Qloc; 626 char filepath[PETSC_MAX_PATH_LEN]; 627 PetscViewer viewer; 628 PetscErrorCode ierr; 629 630 // Set up output 631 PetscFunctionBeginUser; 632 // Print every 'outputfreq' steps 633 if (stepno % user->outputfreq != 0) 634 PetscFunctionReturn(0); 635 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 636 ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 637 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 638 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 639 640 // Output 641 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 642 user->outputdir, stepno + user->contsteps); 643 CHKERRQ(ierr); 644 ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 645 FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 646 ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 647 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 648 if (user->dmviz) { 649 Vec Qrefined, Qrefined_loc; 650 char filepath_refined[PETSC_MAX_PATH_LEN]; 651 PetscViewer viewer_refined; 652 653 ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 654 ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 655 ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 656 CHKERRQ(ierr); 657 ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 658 ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 659 ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 660 CHKERRQ(ierr); 661 ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 662 "%s/nsrefined-%03D.vtu", 663 user->outputdir, stepno + user->contsteps); 664 CHKERRQ(ierr); 665 ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 666 filepath_refined, 667 FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 668 ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 669 ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 670 ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 671 ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 672 } 673 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 674 675 // Save data in a binary file for continuation of simulations 676 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 677 user->outputdir); CHKERRQ(ierr); 678 ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 679 CHKERRQ(ierr); 680 ierr = VecView(Q, viewer); CHKERRQ(ierr); 681 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 682 683 // Save time stamp 684 // Dimensionalize time back 685 time /= user->units->second; 686 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 687 user->outputdir); CHKERRQ(ierr); 688 ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 689 CHKERRQ(ierr); 690 #if PETSC_VERSION_GE(3,13,0) 691 ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 692 #else 693 ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 694 #endif 695 CHKERRQ(ierr); 696 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 697 698 PetscFunctionReturn(0); 699 } 700 701 static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 702 CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 703 CeedElemRestriction restrictq, CeedQFunctionContext ctxSetup, CeedScalar time) { 704 PetscErrorCode ierr; 705 CeedVector multlvec; 706 Vec Multiplicity, MultiplicityLoc; 707 708 SetupContext ctxSetupData; 709 CeedQFunctionContextGetData(ctxSetup, CEED_MEM_HOST, (void **)&ctxSetupData); 710 ctxSetupData->time = time; 711 CeedQFunctionContextRestoreData(ctxSetup, (void **)&ctxSetupData); 712 713 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 714 ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 715 CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 716 ierr = VecZeroEntries(Q); CHKERRQ(ierr); 717 ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 718 719 // Fix multiplicity for output of ICs 720 ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 721 CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 722 ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 723 CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 724 CeedVectorDestroy(&multlvec); 725 ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 726 ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 727 ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 728 CHKERRQ(ierr); 729 ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 730 ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 731 ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 732 ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 733 734 PetscFunctionReturn(0); 735 } 736 737 static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 738 CeedElemRestriction restrictq, CeedBasis basisq, 739 CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 740 PetscErrorCode ierr; 741 CeedQFunction qf_mass; 742 CeedOperator op_mass; 743 CeedVector mceed; 744 Vec Mloc; 745 CeedInt ncompq, qdatasize; 746 747 PetscFunctionBeginUser; 748 CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 749 CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 750 // Create the Q-function that defines the action of the mass operator 751 CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 752 CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 753 CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 754 CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 755 756 // Create the mass operator 757 CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 758 CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 759 CeedOperatorSetField(op_mass, "qdata", restrictqdi, 760 CEED_BASIS_COLLOCATED, qdata); 761 CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 762 763 ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 764 ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 765 CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 766 ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 767 768 { 769 // Compute a lumped mass matrix 770 CeedVector onesvec; 771 CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 772 CeedVectorSetValue(onesvec, 1.0); 773 CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 774 CeedVectorDestroy(&onesvec); 775 CeedOperatorDestroy(&op_mass); 776 CeedVectorDestroy(&mceed); 777 } 778 CeedQFunctionDestroy(&qf_mass); 779 780 ierr = VecZeroEntries(M); CHKERRQ(ierr); 781 ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 782 ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 783 784 // Invert diagonally lumped mass vector for RHS function 785 ierr = VecReciprocal(M); CHKERRQ(ierr); 786 PetscFunctionReturn(0); 787 } 788 789 static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 790 SimpleBC bc, void *ctxSetupData, void *ctxMMS) { 791 PetscErrorCode ierr; 792 793 PetscFunctionBeginUser; 794 { 795 // Configure the finite element space and boundary conditions 796 PetscFE fe; 797 PetscInt ncompq = 5; 798 ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 799 PETSC_FALSE, degree, PETSC_DECIDE, 800 &fe); CHKERRQ(ierr); 801 ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 802 ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 803 ierr = DMCreateDS(dm); CHKERRQ(ierr); 804 if (problem->bc == Exact_Euler) 805 bc->nwall = bc->nslip[0] = bc->nslip[1] = 0; 806 { 807 PetscInt comps[1] = {1}; 808 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 809 1, comps, (void(*)(void))NULL, NULL, bc->nslip[0], 810 bc->slips[0], ctxSetupData); CHKERRQ(ierr); 811 comps[0] = 2; 812 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 813 1, comps, (void(*)(void))NULL, NULL, bc->nslip[1], 814 bc->slips[1], ctxSetupData); CHKERRQ(ierr); 815 comps[0] = 3; 816 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 817 1, comps, (void(*)(void))NULL, NULL, bc->nslip[2], 818 bc->slips[2], ctxSetupData); CHKERRQ(ierr); 819 } 820 if (bc->userbc == PETSC_TRUE) { 821 for (PetscInt c = 0; c < 3; c++) { 822 for (PetscInt s = 0; s < bc->nslip[c]; s++) { 823 for (PetscInt w = 0; w < bc->nwall; w++) { 824 if (bc->slips[c][s] == bc->walls[w]) 825 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 826 "Boundary condition already set on face %D!\n", 827 bc->walls[w]); 828 } 829 } 830 } 831 } 832 // Wall boundary conditions are zero energy density and zero flux for 833 // velocity in advection/advection2d, and zero velocity and zero flux 834 // for mass density and energy density in density_current 835 { 836 if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 837 PetscInt comps[1] = {4}; 838 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 839 1, comps, (void(*)(void))problem->bc, NULL, 840 bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 841 } else if (problem->bc == Exact_DC) { 842 PetscInt comps[3] = {1, 2, 3}; 843 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 844 3, comps, (void(*)(void))problem->bc, NULL, 845 bc->nwall, bc->walls, ctxSetupData); CHKERRQ(ierr); 846 } else if (problem->bc == Exact_Euler) { 847 PetscInt bcMMS[2] = {3, 4}; 848 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", "Face Sets", 0, 849 0, NULL, (void(*)(void))problem->bc, NULL, 850 2, bcMMS, ctxMMS); CHKERRQ(ierr); 851 } else 852 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, 853 "Undefined boundary conditions for this problem"); 854 } 855 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 856 CHKERRQ(ierr); 857 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 858 } 859 { 860 // Empty name for conserved field (because there is only one field) 861 PetscSection section; 862 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 863 ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 864 ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 865 CHKERRQ(ierr); 866 ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 867 CHKERRQ(ierr); 868 ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 869 CHKERRQ(ierr); 870 ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 871 CHKERRQ(ierr); 872 ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 873 CHKERRQ(ierr); 874 } 875 PetscFunctionReturn(0); 876 } 877 878 int main(int argc, char **argv) { 879 PetscInt ierr; 880 MPI_Comm comm; 881 DM dm, dmcoord, dmviz; 882 Mat interpviz; 883 TS ts; 884 TSAdapt adapt; 885 User user; 886 Units units; 887 EulerContext ctxEulerData; 888 char ceedresource[4096] = "/cpu/self"; 889 PetscInt localNelemVol, lnodes, gnodes, steps; 890 const PetscInt ncompq = 5; 891 PetscMPIInt rank; 892 PetscScalar ftime; 893 Vec Q, Qloc, Xloc; 894 Ceed ceed; 895 CeedInt numP, numQ; 896 CeedVector xcorners, qdata, q0ceed; 897 CeedBasis basisx, basisxc, basisq; 898 CeedElemRestriction restrictx, restrictq, restrictqdi; 899 CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol; 900 CeedQFunctionContext ctxSetup, ctxNS, ctxAdvection2d, ctxSurface, ctxEuler; 901 CeedOperator op_setupVol, op_ics; 902 CeedScalar Rd; 903 CeedMemType memtyperequested; 904 PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 905 kgpersquaredms, Joulepercubicm, Joule; 906 problemType problemChoice; 907 problemData *problem = NULL; 908 WindType wind_type; 909 StabilizationType stab; 910 PetscBool implicit; 911 PetscInt viz_refine = 0; 912 struct SimpleBC_ bc = { 913 .nslip = {2, 2, 2}, 914 .slips = {{5, 6}, {3, 4}, {1, 2}} 915 }; 916 double start, cpu_time_used; 917 // Test variables 918 PetscBool test; 919 PetscScalar testtol = 0.; 920 char filepath[PETSC_MAX_PATH_LEN]; 921 // Check PETSc CUDA support 922 PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 923 // *INDENT-OFF* 924 #ifdef PETSC_HAVE_CUDA 925 petschavecuda = PETSC_TRUE; 926 #else 927 petschavecuda = PETSC_FALSE; 928 #endif 929 // *INDENT-ON* 930 931 // Create the libCEED contexts 932 PetscScalar meter = 1e-2; // 1 meter in scaled length units 933 PetscScalar second = 1e-2; // 1 second in scaled time units 934 PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 935 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 936 CeedScalar theta0 = 300.; // K 937 CeedScalar thetaC = -15.; // K 938 CeedScalar P0 = 1.e5; // Pa 939 CeedScalar E_wind = 1.e6; // J 940 CeedScalar N = 0.01; // 1/s 941 CeedScalar cv = 717.; // J/(kg K) 942 CeedScalar cp = 1004.; // J/(kg K) 943 CeedScalar vortex_strength = 5.; // - 944 CeedScalar g = 9.81; // m/s^2 945 CeedScalar lambda = -2./3.; // - 946 CeedScalar mu = 75.; // Pa s, dynamic viscosity 947 // mu = 75 is not physical for air, but is good for numerical stability 948 CeedScalar k = 0.02638; // W/(m K) 949 CeedScalar CtauS = 0.; // dimensionless 950 CeedScalar strong_form = 0.; // [0,1] 951 PetscScalar lx = 8000.; // m 952 PetscScalar ly = 8000.; // m 953 PetscScalar lz = 4000.; // m 954 CeedScalar rc = 1000.; // m (Radius of bubble) 955 PetscScalar resx = 1000.; // m (resolution in x) 956 PetscScalar resy = 1000.; // m (resolution in y) 957 PetscScalar resz = 1000.; // m (resolution in z) 958 PetscInt outputfreq = 10; // - 959 PetscInt contsteps = 0; // - 960 PetscInt degree = 1; // - 961 PetscInt qextra = 2; // - 962 PetscInt qextraSur = 2; // - 963 PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0}; 964 965 ierr = PetscInitialize(&argc, &argv, NULL, help); 966 if (ierr) return ierr; 967 968 // Allocate PETSc context 969 ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 970 ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 971 ierr = PetscMalloc1(1, &ctxEulerData); CHKERRQ(ierr); 972 973 // Parse command line options 974 comm = PETSC_COMM_WORLD; 975 ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 976 NULL); CHKERRQ(ierr); 977 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 978 NULL, ceedresource, ceedresource, 979 sizeof(ceedresource), NULL); CHKERRQ(ierr); 980 ierr = PetscOptionsBool("-test", "Run in test mode", 981 NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr); 982 ierr = PetscOptionsScalar("-compare_final_state_atol", 983 "Test absolute tolerance", 984 NULL, testtol, &testtol, NULL); CHKERRQ(ierr); 985 ierr = PetscOptionsString("-compare_final_state_filename", "Test filename", 986 NULL, filepath, filepath, 987 sizeof(filepath), NULL); CHKERRQ(ierr); 988 problemChoice = NS_DENSITY_CURRENT; 989 ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 990 problemTypes, (PetscEnum)problemChoice, 991 (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 992 problem = &problemOptions[problemChoice]; 993 ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection", 994 NULL, WindTypes, 995 (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION), 996 (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr); 997 PetscInt n = problem->dim; 998 PetscBool userWind; 999 ierr = PetscOptionsRealArray("-problem_advection_wind_translation", 1000 "Constant wind vector", 1001 NULL, wind, &n, &userWind); CHKERRQ(ierr); 1002 if (wind_type == ADVECTION_WIND_ROTATION && userWind) { 1003 ierr = PetscPrintf(comm, 1004 "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n"); 1005 CHKERRQ(ierr); 1006 } 1007 if (wind_type == ADVECTION_WIND_TRANSLATION 1008 && (problemChoice == NS_DENSITY_CURRENT || 1009 problemChoice == NS_EULER_VORTEX)) { 1010 SETERRQ(comm, PETSC_ERR_ARG_INCOMP, 1011 "-problem_advection_wind translation is not defined for -problem density_current or -problem euler_vortex"); 1012 } 1013 ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 1014 StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 1015 (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 1016 ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 1017 NULL, implicit=PETSC_FALSE, &implicit, NULL); 1018 CHKERRQ(ierr); 1019 if (!implicit && stab != STAB_NONE) { 1020 ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 1021 CHKERRQ(ierr); 1022 } 1023 { 1024 PetscInt len; 1025 PetscBool flg; 1026 ierr = PetscOptionsIntArray("-bc_wall", 1027 "Use wall boundary conditions on this list of faces", 1028 NULL, bc.walls, 1029 (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 1030 &len), &flg); CHKERRQ(ierr); 1031 if (flg) { 1032 bc.nwall = len; 1033 // Using a no-slip wall disables automatic slip walls (they must be set explicitly) 1034 bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0; 1035 } 1036 for (PetscInt j=0; j<3; j++) { 1037 const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1038 ierr = PetscOptionsIntArray(flags[j], 1039 "Use slip boundary conditions on this list of faces", 1040 NULL, bc.slips[j], 1041 (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 1042 &len), &flg); 1043 CHKERRQ(ierr); 1044 if (flg) { 1045 bc.nslip[j] = len; 1046 bc.userbc = PETSC_TRUE; 1047 } 1048 } 1049 } 1050 ierr = PetscOptionsInt("-viz_refine", 1051 "Regular refinement levels for visualization", 1052 NULL, viz_refine, &viz_refine, NULL); 1053 CHKERRQ(ierr); 1054 ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 1055 NULL, meter, &meter, NULL); CHKERRQ(ierr); 1056 meter = fabs(meter); 1057 ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 1058 NULL, second, &second, NULL); CHKERRQ(ierr); 1059 second = fabs(second); 1060 ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 1061 NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 1062 kilogram = fabs(kilogram); 1063 ierr = PetscOptionsScalar("-units_Kelvin", 1064 "1 Kelvin in scaled temperature units", 1065 NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 1066 Kelvin = fabs(Kelvin); 1067 ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 1068 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 1069 ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 1070 NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 1071 ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 1072 NULL, P0, &P0, NULL); CHKERRQ(ierr); 1073 ierr = PetscOptionsScalar("-E_wind", "Total energy of inflow wind", 1074 NULL, E_wind, &E_wind, NULL); CHKERRQ(ierr); 1075 ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 1076 NULL, N, &N, NULL); CHKERRQ(ierr); 1077 ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 1078 NULL, cv, &cv, NULL); CHKERRQ(ierr); 1079 ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 1080 NULL, cp, &cp, NULL); CHKERRQ(ierr); 1081 PetscBool userVortex; 1082 ierr = PetscOptionsScalar("-vortex_strength", "Strength of Vortex", 1083 NULL, vortex_strength, &vortex_strength, &userVortex); 1084 CHKERRQ(ierr); 1085 if (problemChoice != NS_EULER_VORTEX && userVortex) { 1086 ierr = PetscPrintf(comm, 1087 "Warning! Use -vortex_strength only with -problem euler_vortex\n"); 1088 CHKERRQ(ierr); 1089 } 1090 ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 1091 NULL, g, &g, NULL); CHKERRQ(ierr); 1092 ierr = PetscOptionsScalar("-lambda", 1093 "Stokes hypothesis second viscosity coefficient", 1094 NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 1095 ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 1096 NULL, mu, &mu, NULL); CHKERRQ(ierr); 1097 ierr = PetscOptionsScalar("-k", "Thermal conductivity", 1098 NULL, k, &k, NULL); CHKERRQ(ierr); 1099 ierr = PetscOptionsScalar("-CtauS", 1100 "Scale coefficient for tau (nondimensional)", 1101 NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 1102 if (stab == STAB_NONE && CtauS != 0) { 1103 ierr = PetscPrintf(comm, 1104 "Warning! Use -CtauS only with -stab su or -stab supg\n"); 1105 CHKERRQ(ierr); 1106 } 1107 ierr = PetscOptionsScalar("-strong_form", 1108 "Strong (1) or weak/integrated by parts (0) advection residual", 1109 NULL, strong_form, &strong_form, NULL); 1110 CHKERRQ(ierr); 1111 if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) { 1112 ierr = PetscPrintf(comm, 1113 "Warning! Problem density_current does not support -CtauS or -strong_form\n"); 1114 CHKERRQ(ierr); 1115 } 1116 ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 1117 NULL, lx, &lx, NULL); CHKERRQ(ierr); 1118 ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 1119 NULL, ly, &ly, NULL); CHKERRQ(ierr); 1120 ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 1121 NULL, lz, &lz, NULL); CHKERRQ(ierr); 1122 ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 1123 NULL, rc, &rc, NULL); CHKERRQ(ierr); 1124 ierr = PetscOptionsScalar("-resx","Target resolution in x", 1125 NULL, resx, &resx, NULL); CHKERRQ(ierr); 1126 ierr = PetscOptionsScalar("-resy","Target resolution in y", 1127 NULL, resy, &resy, NULL); CHKERRQ(ierr); 1128 ierr = PetscOptionsScalar("-resz","Target resolution in z", 1129 NULL, resz, &resz, NULL); CHKERRQ(ierr); 1130 n = problem->dim; 1131 center[0] = 0.5 * lx; 1132 center[1] = 0.5 * ly; 1133 center[2] = 0.5 * lz; 1134 ierr = PetscOptionsRealArray("-center", "Location of bubble center", 1135 NULL, center, &n, NULL); CHKERRQ(ierr); 1136 n = problem->dim; 1137 ierr = PetscOptionsRealArray("-dc_axis", 1138 "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 1139 NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 1140 { 1141 PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 1142 PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 1143 if (norm > 0) { 1144 for (int i=0; i<3; i++) dc_axis[i] /= norm; 1145 } 1146 } 1147 ierr = PetscOptionsInt("-output_freq", 1148 "Frequency of output, in number of steps", 1149 NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 1150 ierr = PetscOptionsInt("-continue", "Continue from previous solution", 1151 NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 1152 ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 1153 NULL, degree, °ree, NULL); CHKERRQ(ierr); 1154 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 1155 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 1156 PetscBool userQextraSur; 1157 ierr = PetscOptionsInt("-qextra_boundary", 1158 "Number of extra quadrature points on in/outflow faces", 1159 NULL, qextraSur, &qextraSur, &userQextraSur); 1160 CHKERRQ(ierr); 1161 if ((wind_type == ADVECTION_WIND_ROTATION 1162 || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) { 1163 ierr = PetscPrintf(comm, 1164 "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n"); 1165 CHKERRQ(ierr); 1166 } 1167 ierr = PetscStrncpy(user->outputdir, ".", 2); CHKERRQ(ierr); 1168 ierr = PetscOptionsString("-output_dir", "Output directory", 1169 NULL, user->outputdir, user->outputdir, 1170 sizeof(user->outputdir), NULL); CHKERRQ(ierr); 1171 memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 1172 ierr = PetscOptionsEnum("-memtype", 1173 "CEED MemType requested", NULL, 1174 memTypes, (PetscEnum)memtyperequested, 1175 (PetscEnum *)&memtyperequested, &setmemtyperequest); 1176 CHKERRQ(ierr); 1177 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 1178 1179 // Define derived units 1180 Pascal = kilogram / (meter * PetscSqr(second)); 1181 JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 1182 mpersquareds = meter / PetscSqr(second); 1183 WpermK = kilogram * meter / (pow(second,3) * Kelvin); 1184 kgpercubicm = kilogram / pow(meter,3); 1185 kgpersquaredms = kilogram / (PetscSqr(meter) * second); 1186 Joulepercubicm = kilogram / (meter * PetscSqr(second)); 1187 Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 1188 1189 // Scale variables to desired units 1190 theta0 *= Kelvin; 1191 thetaC *= Kelvin; 1192 P0 *= Pascal; 1193 E_wind *= Joule; 1194 N *= (1./second); 1195 cv *= JperkgK; 1196 cp *= JperkgK; 1197 Rd = cp - cv; 1198 g *= mpersquareds; 1199 mu *= Pascal * second; 1200 k *= WpermK; 1201 lx = fabs(lx) * meter; 1202 ly = fabs(ly) * meter; 1203 lz = fabs(lz) * meter; 1204 rc = fabs(rc) * meter; 1205 resx = fabs(resx) * meter; 1206 resy = fabs(resy) * meter; 1207 resz = fabs(resz) * meter; 1208 for (int i=0; i<3; i++) center[i] *= meter; 1209 1210 const CeedInt dim = problem->dim, ncompx = problem->dim, 1211 qdatasizeVol = problem->qdatasizeVol; 1212 // Set up the libCEED context 1213 struct SetupContext_ ctxSetupData = { 1214 .theta0 = theta0, 1215 .thetaC = thetaC, 1216 .P0 = P0, 1217 .N = N, 1218 .cv = cv, 1219 .cp = cp, 1220 .Rd = Rd, 1221 .g = g, 1222 .rc = rc, 1223 .lx = lx, 1224 .ly = ly, 1225 .lz = lz, 1226 .center[0] = center[0], 1227 .center[1] = center[1], 1228 .center[2] = center[2], 1229 .dc_axis[0] = dc_axis[0], 1230 .dc_axis[1] = dc_axis[1], 1231 .dc_axis[2] = dc_axis[2], 1232 .wind[0] = wind[0], 1233 .wind[1] = wind[1], 1234 .wind[2] = wind[2], 1235 .time = 0, 1236 .vortex_strength = vortex_strength, 1237 .wind_type = wind_type, 1238 }; 1239 1240 // Create the mesh 1241 { 1242 const PetscReal scale[3] = {lx, ly, lz}; 1243 ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 1244 NULL, PETSC_TRUE, &dm); 1245 CHKERRQ(ierr); 1246 } 1247 1248 // Distribute the mesh over processes 1249 { 1250 DM dmDist = NULL; 1251 PetscPartitioner part; 1252 1253 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 1254 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 1255 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 1256 if (dmDist) { 1257 ierr = DMDestroy(&dm); CHKERRQ(ierr); 1258 dm = dmDist; 1259 } 1260 } 1261 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 1262 1263 // Setup DM 1264 ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 1265 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 1266 ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetupData, ctxEulerData); 1267 CHKERRQ(ierr); 1268 1269 // Refine DM for high-order viz 1270 dmviz = NULL; 1271 interpviz = NULL; 1272 if (viz_refine) { 1273 DM dmhierarchy[viz_refine+1]; 1274 1275 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 1276 dmhierarchy[0] = dm; 1277 for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 1278 Mat interp_next; 1279 1280 ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 1281 CHKERRQ(ierr); 1282 ierr = DMClearDS(dmhierarchy[i+1]); CHKERRQ(ierr); 1283 ierr = DMClearFields(dmhierarchy[i+1]); CHKERRQ(ierr); 1284 ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 1285 d = (d + 1) / 2; 1286 if (i + 1 == viz_refine) d = 1; 1287 ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetupData, 1288 ctxEulerData); CHKERRQ(ierr); 1289 ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 1290 &interp_next, NULL); CHKERRQ(ierr); 1291 if (!i) interpviz = interp_next; 1292 else { 1293 Mat C; 1294 ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 1295 PETSC_DECIDE, &C); CHKERRQ(ierr); 1296 ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 1297 ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1298 interpviz = C; 1299 } 1300 } 1301 for (PetscInt i=1; i<viz_refine; i++) { 1302 ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 1303 } 1304 dmviz = dmhierarchy[viz_refine]; 1305 } 1306 ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 1307 ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 1308 ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 1309 lnodes /= ncompq; 1310 1311 // Initialize CEED 1312 CeedInit(ceedresource, &ceed); 1313 // Set memtype 1314 CeedMemType memtypebackend; 1315 CeedGetPreferredMemType(ceed, &memtypebackend); 1316 // Check memtype compatibility 1317 if (!setmemtyperequest) 1318 memtyperequested = memtypebackend; 1319 else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 1320 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 1321 "PETSc was not built with CUDA. " 1322 "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 1323 1324 // Set number of 1D nodes and quadrature points 1325 numP = degree + 1; 1326 numQ = numP + qextra; 1327 1328 // Print summary 1329 if (!test) { 1330 CeedInt gdofs, odofs; 1331 int comm_size; 1332 char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 1333 ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 1334 ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 1335 gnodes = gdofs/ncompq; 1336 ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 1337 ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 1338 sizeof(box_faces_str), NULL); CHKERRQ(ierr); 1339 const char *usedresource; 1340 CeedGetResource(ceed, &usedresource); 1341 1342 ierr = PetscPrintf(comm, 1343 "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 1344 " rank(s) : %d\n" 1345 " Problem:\n" 1346 " Problem Name : %s\n" 1347 " Stabilization : %s\n" 1348 " PETSc:\n" 1349 " Box Faces : %s\n" 1350 " libCEED:\n" 1351 " libCEED Backend : %s\n" 1352 " libCEED Backend MemType : %s\n" 1353 " libCEED User Requested MemType : %s\n" 1354 " Mesh:\n" 1355 " Number of 1D Basis Nodes (P) : %d\n" 1356 " Number of 1D Quadrature Points (Q) : %d\n" 1357 " Global DoFs : %D\n" 1358 " Owned DoFs : %D\n" 1359 " DoFs per node : %D\n" 1360 " Global nodes : %D\n" 1361 " Owned nodes : %D\n", 1362 comm_size, problemTypes[problemChoice], 1363 StabilizationTypes[stab], box_faces_str, usedresource, 1364 CeedMemTypes[memtypebackend], 1365 (setmemtyperequest) ? 1366 CeedMemTypes[memtyperequested] : "none", 1367 numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes); 1368 CHKERRQ(ierr); 1369 } 1370 1371 // Set up global mass vector 1372 ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 1373 1374 // Set up libCEED 1375 // CEED Bases 1376 CeedInit(ceedresource, &ceed); 1377 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 1378 &basisq); 1379 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 1380 &basisx); 1381 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 1382 CEED_GAUSS_LOBATTO, &basisxc); 1383 ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 1384 ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 1385 CHKERRQ(ierr); 1386 1387 // CEED Restrictions 1388 ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ, 1389 qdatasizeVol, &restrictq, &restrictx, 1390 &restrictqdi); CHKERRQ(ierr); 1391 1392 ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1393 ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1394 1395 // Create the CEED vectors that will be needed in setup 1396 CeedInt NqptsVol; 1397 CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol); 1398 CeedElemRestrictionGetNumElements(restrictq, &localNelemVol); 1399 CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata); 1400 CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1401 1402 // Create the Q-function that builds the quadrature data for the NS operator 1403 CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc, 1404 &qf_setupVol); 1405 CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD); 1406 CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT); 1407 CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1408 1409 // Create the Q-function that sets the ICs of the operator 1410 CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1411 CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1412 CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1413 1414 qf_rhsVol = NULL; 1415 if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator 1416 CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs, 1417 problem->applyVol_rhs_loc, &qf_rhsVol); 1418 CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP); 1419 CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 1420 CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1421 CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP); 1422 CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP); 1423 CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1424 } 1425 1426 qf_ifunctionVol = NULL; 1427 if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction 1428 CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction, 1429 problem->applyVol_ifunction_loc, &qf_ifunctionVol); 1430 CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP); 1431 CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 1432 CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP); 1433 CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1434 CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP); 1435 CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP); 1436 CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1437 } 1438 1439 // Create the operator that builds the quadrature data for the NS operator 1440 CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol); 1441 CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1442 CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE, 1443 basisx, CEED_VECTOR_NONE); 1444 CeedOperatorSetField(op_setupVol, "qdata", restrictqdi, 1445 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1446 1447 // Create the operator that sets the ICs 1448 CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 1449 CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 1450 CeedOperatorSetField(op_ics, "q0", restrictq, 1451 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1452 1453 CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 1454 CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 1455 CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1456 1457 if (qf_rhsVol) { // Create the RHS physics operator 1458 CeedOperator op; 1459 CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op); 1460 CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1461 CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1462 CeedOperatorSetField(op, "qdata", restrictqdi, 1463 CEED_BASIS_COLLOCATED, qdata); 1464 CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1465 CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1466 CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1467 user->op_rhs_vol = op; 1468 } 1469 1470 if (qf_ifunctionVol) { // Create the IFunction operator 1471 CeedOperator op; 1472 CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op); 1473 CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1474 CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1475 CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 1476 CeedOperatorSetField(op, "qdata", restrictqdi, 1477 CEED_BASIS_COLLOCATED, qdata); 1478 CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1479 CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1480 CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1481 user->op_ifunction_vol = op; 1482 } 1483 1484 // Set up CEED for the boundaries 1485 CeedInt height = 1; 1486 CeedInt dimSur = dim - height; 1487 CeedInt numP_Sur = degree + 1; 1488 CeedInt numQ_Sur = numP_Sur + qextraSur; 1489 const CeedInt qdatasizeSur = problem->qdatasizeSur; 1490 CeedBasis basisxSur, basisxcSur, basisqSur; 1491 CeedInt NqptsSur; 1492 CeedQFunction qf_setupSur, qf_applySur; 1493 1494 // CEED bases for the boundaries 1495 CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, 1496 CEED_GAUSS, 1497 &basisqSur); 1498 CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS, 1499 &basisxSur); 1500 CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur, 1501 CEED_GAUSS_LOBATTO, &basisxcSur); 1502 CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur); 1503 1504 // Create the Q-function that builds the quadrature data for the Surface operator 1505 CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc, 1506 &qf_setupSur); 1507 CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD); 1508 CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT); 1509 CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 1510 1511 // Creat Q-Function for Boundaries 1512 // -- Defined for Advection(2d) test cases 1513 qf_applySur = NULL; 1514 if (problem->applySur) { 1515 CeedQFunctionCreateInterior(ceed, 1, problem->applySur, 1516 problem->applySur_loc, &qf_applySur); 1517 CeedQFunctionAddInput(qf_applySur, "q", ncompq, CEED_EVAL_INTERP); 1518 CeedQFunctionAddInput(qf_applySur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 1519 CeedQFunctionAddInput(qf_applySur, "x", ncompx, CEED_EVAL_INTERP); 1520 CeedQFunctionAddOutput(qf_applySur, "v", ncompq, CEED_EVAL_INTERP); 1521 } 1522 1523 // Create CEED Operator for the whole domain 1524 if (!implicit) 1525 ierr = CreateOperatorForDomain(ceed, dm, &bc, problemChoice, wind_type, 1526 user->op_rhs_vol, qf_applySur, 1527 qf_setupSur, height, numP_Sur, numQ_Sur, 1528 qdatasizeSur, NqptsSur, basisxSur, 1529 basisqSur, &user->op_rhs); 1530 CHKERRQ(ierr); 1531 if (implicit) 1532 ierr = CreateOperatorForDomain(ceed, dm, &bc, problemChoice, wind_type, 1533 user->op_ifunction_vol, qf_applySur, 1534 qf_setupSur, height, numP_Sur, numQ_Sur, 1535 qdatasizeSur, NqptsSur, basisxSur, 1536 basisqSur, &user->op_ifunction); 1537 CHKERRQ(ierr); 1538 // Set up contex for QFunctions 1539 CeedQFunctionContextCreate(ceed, &ctxSetup); 1540 CeedQFunctionContextSetData(ctxSetup, CEED_MEM_HOST, CEED_USE_POINTER, 1541 sizeof ctxSetupData, &ctxSetupData); 1542 if (qf_ics && problemChoice != NS_EULER_VORTEX) 1543 CeedQFunctionSetContext(qf_ics, ctxSetup); 1544 1545 CeedScalar ctxNSData[8] = {lambda, mu, k, cv, cp, g, Rd}; 1546 CeedQFunctionContextCreate(ceed, &ctxNS); 1547 CeedQFunctionContextSetData(ctxNS, CEED_MEM_HOST, CEED_USE_POINTER, 1548 sizeof ctxNSData, &ctxNSData); 1549 1550 struct Advection2dContext_ ctxAdvection2d = { 1551 .CtauS = CtauS, 1552 .strong_form = strong_form, 1553 .stabilization = stab, 1554 }; 1555 CeedQFunctionContextCreate(ceed, &ctxAdvection2d); 1556 CeedQFunctionContextSetData(ctxAdvection2d, CEED_MEM_HOST, CEED_USE_POINTER, 1557 sizeof ctxAdvection2dData, &ctxAdvection2dData); 1558 1559 struct SurfaceContext_ ctxSurfaceData = { 1560 .E_wind = E_wind, 1561 .strong_form = strong_form, 1562 .implicit = implicit, 1563 }; 1564 CeedQFunctionContextCreate(ceed, &ctxSurface); 1565 CeedQFunctionContextSetData(ctxSurface, CEED_MEM_HOST, CEED_USE_POINTER, 1566 sizeof ctxSurfaceData, &ctxSurfaceData); 1567 1568 // Set up ctxEulerData structure 1569 ctxEulerData->time = 0.; 1570 ctxEulerData->currentTime = 0.; 1571 ctxEulerData->wind[0] = wind[0]; 1572 ctxEulerData->wind[1] = wind[1]; 1573 ctxEulerData->wind[2] = wind[2]; 1574 ctxEulerData->center[0] = center[0]; 1575 ctxEulerData->center[1] = center[1]; 1576 ctxEulerData->center[2] = center[2]; 1577 ctxEulerData->vortex_strength = vortex_strength; 1578 1579 user->ctxEulerData = ctxEulerData; 1580 1581 CeedQFunctionContextCreate(ceed, &ctxEuler); 1582 CeedQFunctionContextSetData(ctxEuler, CEED_MEM_HOST, CEED_USE_POINTER, 1583 sizeof *ctxEulerData, ctxEulerData); // ToDo: check the pointer 1584 1585 switch (problemChoice) { 1586 case NS_DENSITY_CURRENT: 1587 if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, ctxNS); 1588 if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, ctxNS); 1589 break; 1590 case NS_ADVECTION: 1591 case NS_ADVECTION2D: 1592 if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, ctxAdvection2d); 1593 if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, ctxAdvection2d); 1594 if (qf_applySur) CeedQFunctionSetContext(qf_applySur, ctxSurface); 1595 case NS_EULER_VORTEX: 1596 if (qf_ics) CeedQFunctionSetContext(qf_ics, ctxEuler); // ToDo: check the pointer 1597 if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, ctxEuler); 1598 if (qf_applySur) CeedQFunctionSetContext(qf_applySur, ctxEuler); 1599 } 1600 1601 // Set up PETSc context 1602 // Set up units structure 1603 units->meter = meter; 1604 units->kilogram = kilogram; 1605 units->second = second; 1606 units->Kelvin = Kelvin; 1607 units->Pascal = Pascal; 1608 units->JperkgK = JperkgK; 1609 units->mpersquareds = mpersquareds; 1610 units->WpermK = WpermK; 1611 units->kgpercubicm = kgpercubicm; 1612 units->kgpersquaredms = kgpersquaredms; 1613 units->Joulepercubicm = Joulepercubicm; 1614 units->Joule = Joule; 1615 1616 // Set up user structure 1617 user->comm = comm; 1618 user->outputfreq = outputfreq; 1619 user->contsteps = contsteps; 1620 user->units = units; 1621 user->dm = dm; 1622 user->dmviz = dmviz; 1623 user->interpviz = interpviz; 1624 user->ceed = ceed; 1625 1626 // Calculate qdata and ICs 1627 // Set up state global and local vectors 1628 ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1629 1630 ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1631 1632 // Apply Setup Ceed Operators 1633 ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 1634 CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 1635 ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1636 user->M); CHKERRQ(ierr); 1637 1638 ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 1639 ctxSetup, 0.0); CHKERRQ(ierr); 1640 if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1641 // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1642 // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1643 // slower execution. 1644 Vec Qbc; 1645 ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1646 ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1647 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1648 ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1649 ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1650 ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1651 ierr = PetscObjectComposeFunction((PetscObject)dm, 1652 "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 1653 CHKERRQ(ierr); 1654 } 1655 1656 MPI_Comm_rank(comm, &rank); 1657 if (!rank) {ierr = PetscMkdir(user->outputdir); CHKERRQ(ierr);} 1658 // Gather initial Q values 1659 // In case of continuation of simulation, set up initial values from binary file 1660 if (contsteps) { // continue from existent solution 1661 PetscViewer viewer; 1662 char filepath[PETSC_MAX_PATH_LEN]; 1663 // Read input 1664 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1665 user->outputdir); 1666 CHKERRQ(ierr); 1667 ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1668 CHKERRQ(ierr); 1669 ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1670 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1671 } 1672 ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1673 1674 // Create and setup TS 1675 ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1676 ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1677 if (implicit) { 1678 ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1679 if (user->op_ifunction) { 1680 ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1681 } else { // Implicit integrators can fall back to using an RHSFunction 1682 ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1683 } 1684 } else { 1685 if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, 1686 "Problem does not provide RHSFunction"); 1687 ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1688 ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1689 ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1690 } 1691 ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1692 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1693 ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 1694 if (test) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);} 1695 ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1696 ierr = TSAdaptSetStepLimits(adapt, 1697 1.e-12 * units->second, 1698 1.e2 * units->second); CHKERRQ(ierr); 1699 ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1700 if (!contsteps) { // print initial condition 1701 if (!test) { 1702 ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1703 } 1704 } else { // continue from time of last output 1705 PetscReal time; 1706 PetscInt count; 1707 PetscViewer viewer; 1708 char filepath[PETSC_MAX_PATH_LEN]; 1709 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1710 user->outputdir); CHKERRQ(ierr); 1711 ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1712 CHKERRQ(ierr); 1713 ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1714 CHKERRQ(ierr); 1715 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1716 ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1717 } 1718 if (!test) { 1719 ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1720 } 1721 1722 // Solve 1723 start = MPI_Wtime(); 1724 ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1725 ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1726 cpu_time_used = MPI_Wtime() - start; 1727 ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1728 ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1729 comm); CHKERRQ(ierr); 1730 if (!test) { 1731 ierr = PetscPrintf(PETSC_COMM_WORLD, 1732 "Time taken for solution (sec): %g\n", 1733 (double)cpu_time_used); CHKERRQ(ierr); 1734 } 1735 1736 // Get error 1737 if (problem->non_zero_time && !test) { 1738 Vec Qexact, Qexactloc; 1739 PetscReal norm; 1740 ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1741 ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1742 ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1743 1744 ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 1745 restrictq, ctxSetup, ftime); CHKERRQ(ierr); 1746 1747 ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1748 ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1749 CeedVectorDestroy(&q0ceed); 1750 ierr = PetscPrintf(PETSC_COMM_WORLD, 1751 "Max Error: %g\n", 1752 (double)norm); CHKERRQ(ierr); 1753 // Clean up vectors 1754 ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1755 ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1756 } 1757 1758 // Output Statistics 1759 ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr); 1760 if (!test) { 1761 ierr = PetscPrintf(PETSC_COMM_WORLD, 1762 "Time integrator took %D time steps to reach final time %g\n", 1763 steps, (double)ftime); CHKERRQ(ierr); 1764 } 1765 // Output numerical values from command line 1766 ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr); 1767 1768 // Compare reference solution values with current test run for CI 1769 if (test) { 1770 PetscViewer viewer; 1771 // Read reference file 1772 Vec Qref; 1773 PetscReal error, Qrefnorm; 1774 ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr); 1775 ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1776 CHKERRQ(ierr); 1777 ierr = VecLoad(Qref, viewer); CHKERRQ(ierr); 1778 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1779 1780 // Compute error with respect to reference solution 1781 ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr); 1782 ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr); 1783 ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr); 1784 ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr); 1785 ierr = VecDestroy(&Qref); CHKERRQ(ierr); 1786 // Check error 1787 if (error > testtol) { 1788 ierr = PetscPrintf(PETSC_COMM_WORLD, 1789 "Test failed with error norm %g\n", 1790 (double)error); CHKERRQ(ierr); 1791 } 1792 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1793 } 1794 1795 // Clean up libCEED 1796 CeedVectorDestroy(&qdata); 1797 CeedVectorDestroy(&user->qceed); 1798 CeedVectorDestroy(&user->qdotceed); 1799 CeedVectorDestroy(&user->gceed); 1800 CeedVectorDestroy(&xcorners); 1801 CeedBasisDestroy(&basisq); 1802 CeedBasisDestroy(&basisx); 1803 CeedBasisDestroy(&basisxc); 1804 CeedElemRestrictionDestroy(&restrictq); 1805 CeedElemRestrictionDestroy(&restrictx); 1806 CeedElemRestrictionDestroy(&restrictqdi); 1807 CeedQFunctionDestroy(&qf_setupVol); 1808 CeedQFunctionDestroy(&qf_ics); 1809 CeedQFunctionDestroy(&qf_rhsVol); 1810 CeedQFunctionDestroy(&qf_ifunctionVol); 1811 CeedQFunctionContextDestroy(&ctxSetup); 1812 CeedQFunctionContextDestroy(&ctxNS); 1813 CeedQFunctionContextDestroy(&ctxAdvection2d); 1814 CeedQFunctionContextDestroy(&ctxSurface); 1815 CeedOperatorDestroy(&op_setupVol); 1816 CeedOperatorDestroy(&op_ics); 1817 CeedOperatorDestroy(&user->op_rhs_vol); 1818 CeedOperatorDestroy(&user->op_ifunction_vol); 1819 CeedDestroy(&ceed); 1820 CeedBasisDestroy(&basisqSur); 1821 CeedBasisDestroy(&basisxSur); 1822 CeedBasisDestroy(&basisxcSur); 1823 CeedQFunctionDestroy(&qf_setupSur); 1824 CeedQFunctionDestroy(&qf_applySur); 1825 CeedOperatorDestroy(&user->op_rhs); 1826 CeedOperatorDestroy(&user->op_ifunction); 1827 1828 // Clean up PETSc 1829 ierr = VecDestroy(&Q); CHKERRQ(ierr); 1830 ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1831 ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1832 ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1833 ierr = TSDestroy(&ts); CHKERRQ(ierr); 1834 ierr = DMDestroy(&dm); CHKERRQ(ierr); 1835 ierr = PetscFree(units); CHKERRQ(ierr); 1836 ierr = PetscFree(user); CHKERRQ(ierr); 1837 ierr = PetscFree(ctxEulerData); CHKERRQ(ierr); 1838 return PetscFinalize(); 1839 } 1840