1 // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. 2 // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. 3 // All Rights 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 #include "ceed-ref.h" 18 19 //------------------------------------------------------------------------------ 20 // Setup Input/Output Fields 21 //------------------------------------------------------------------------------ 22 static int CeedOperatorSetupFields_Ref(CeedQFunction qf, CeedOperator op, 23 bool inOrOut, 24 CeedVector *fullevecs, CeedVector *evecs, 25 CeedVector *qvecs, CeedInt starte, 26 CeedInt numfields, CeedInt Q) { 27 CeedInt dim, ierr, size, P; 28 Ceed ceed; 29 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 30 CeedBasis basis; 31 CeedElemRestriction Erestrict; 32 CeedOperatorField *opfields; 33 CeedQFunctionField *qffields; 34 if (inOrOut) { 35 ierr = CeedOperatorGetFields(op, NULL, &opfields); CeedChk(ierr); 36 ierr = CeedQFunctionGetFields(qf, NULL, &qffields); CeedChk(ierr); 37 } else { 38 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 39 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 40 } 41 42 // Loop over fields 43 for (CeedInt i=0; i<numfields; i++) { 44 CeedEvalMode emode; 45 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 46 47 if (emode != CEED_EVAL_WEIGHT) { 48 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict); 49 CeedChk(ierr); 50 ierr = CeedElemRestrictionCreateVector(Erestrict, NULL, 51 &fullevecs[i+starte]); 52 CeedChk(ierr); 53 } 54 55 switch(emode) { 56 case CEED_EVAL_NONE: 57 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 58 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 59 break; 60 case CEED_EVAL_INTERP: 61 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 62 ierr = CeedElemRestrictionGetElementSize(Erestrict, &P); 63 CeedChk(ierr); 64 ierr = CeedVectorCreate(ceed, P*size, &evecs[i]); CeedChk(ierr); 65 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 66 break; 67 case CEED_EVAL_GRAD: 68 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 69 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 70 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 71 ierr = CeedElemRestrictionGetElementSize(Erestrict, &P); 72 CeedChk(ierr); 73 ierr = CeedVectorCreate(ceed, P*size/dim, &evecs[i]); CeedChk(ierr); 74 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 75 break; 76 case CEED_EVAL_WEIGHT: // Only on input fields 77 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 78 ierr = CeedVectorCreate(ceed, Q, &qvecs[i]); CeedChk(ierr); 79 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, 80 CEED_VECTOR_NONE, qvecs[i]); CeedChk(ierr); 81 break; 82 case CEED_EVAL_DIV: 83 break; // Not implemented 84 case CEED_EVAL_CURL: 85 break; // Not implemented 86 } 87 } 88 return 0; 89 } 90 91 //------------------------------------------------------------------------------ 92 // Setup Operator 93 //------------------------------------------------------------------------------/* 94 static int CeedOperatorSetup_Ref(CeedOperator op) { 95 int ierr; 96 bool setupdone; 97 ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChk(ierr); 98 if (setupdone) return 0; 99 Ceed ceed; 100 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 101 CeedOperator_Ref *impl; 102 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 103 CeedQFunction qf; 104 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 105 CeedInt Q, numinputfields, numoutputfields; 106 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 107 ierr = CeedQFunctionIsIdentity(qf, &impl->identityqf); CeedChk(ierr); 108 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 109 CeedChk(ierr); 110 CeedOperatorField *opinputfields, *opoutputfields; 111 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 112 CeedChk(ierr); 113 CeedQFunctionField *qfinputfields, *qfoutputfields; 114 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 115 CeedChk(ierr); 116 117 // Allocate 118 ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs); 119 CeedChk(ierr); 120 ierr = CeedCalloc(numinputfields + numoutputfields, &impl->edata); 121 CeedChk(ierr); 122 123 ierr = CeedCalloc(16, &impl->inputstate); CeedChk(ierr); 124 ierr = CeedCalloc(16, &impl->evecsin); CeedChk(ierr); 125 ierr = CeedCalloc(16, &impl->evecsout); CeedChk(ierr); 126 ierr = CeedCalloc(16, &impl->qvecsin); CeedChk(ierr); 127 ierr = CeedCalloc(16, &impl->qvecsout); CeedChk(ierr); 128 129 impl->numein = numinputfields; impl->numeout = numoutputfields; 130 131 // Set up infield and outfield evecs and qvecs 132 // Infields 133 ierr = CeedOperatorSetupFields_Ref(qf, op, 0, impl->evecs, 134 impl->evecsin, impl->qvecsin, 0, 135 numinputfields, Q); 136 CeedChk(ierr); 137 // Outfields 138 ierr = CeedOperatorSetupFields_Ref(qf, op, 1, impl->evecs, 139 impl->evecsout, impl->qvecsout, 140 numinputfields, numoutputfields, Q); 141 CeedChk(ierr); 142 143 // Identity QFunctions 144 if (impl->identityqf) { 145 CeedEvalMode inmode, outmode; 146 CeedQFunctionField *infields, *outfields; 147 ierr = CeedQFunctionGetFields(qf, &infields, &outfields); CeedChk(ierr); 148 149 for (CeedInt i=0; i<numinputfields; i++) { 150 ierr = CeedQFunctionFieldGetEvalMode(infields[i], &inmode); 151 CeedChk(ierr); 152 ierr = CeedQFunctionFieldGetEvalMode(outfields[i], &outmode); 153 CeedChk(ierr); 154 155 ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr); 156 impl->qvecsout[i] = impl->qvecsin[i]; 157 ierr = CeedVectorAddReference(impl->qvecsin[i]); CeedChk(ierr); 158 } 159 } 160 161 ierr = CeedOperatorSetSetupDone(op); CeedChk(ierr); 162 163 return 0; 164 } 165 166 //------------------------------------------------------------------------------ 167 // Setup Operator Inputs 168 //------------------------------------------------------------------------------ 169 static inline int CeedOperatorSetupInputs_Ref(CeedInt numinputfields, 170 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 171 CeedVector invec, const bool skipactive, CeedOperator_Ref *impl, 172 CeedRequest *request) { 173 CeedInt ierr; 174 CeedEvalMode emode; 175 CeedVector vec; 176 CeedElemRestriction Erestrict; 177 uint64_t state; 178 179 for (CeedInt i=0; i<numinputfields; i++) { 180 // Get input vector 181 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 182 if (vec == CEED_VECTOR_ACTIVE) { 183 if (skipactive) 184 continue; 185 else 186 vec = invec; 187 } 188 189 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 190 CeedChk(ierr); 191 // Restrict and Evec 192 if (emode == CEED_EVAL_WEIGHT) { // Skip 193 } else { 194 // Restrict 195 ierr = CeedVectorGetState(vec, &state); CeedChk(ierr); 196 // Skip restriction if input is unchanged 197 if (state != impl->inputstate[i] || vec == invec) { 198 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 199 CeedChk(ierr); 200 ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec, 201 impl->evecs[i], request); CeedChk(ierr); 202 impl->inputstate[i] = state; 203 } 204 // Get evec 205 ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_HOST, 206 (const CeedScalar **) &impl->edata[i]); 207 CeedChk(ierr); 208 } 209 } 210 return 0; 211 } 212 213 //------------------------------------------------------------------------------ 214 // Input Basis Action 215 //------------------------------------------------------------------------------ 216 static inline int CeedOperatorInputBasis_Ref(CeedInt e, CeedInt Q, 217 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 218 CeedInt numinputfields, const bool skipactive, CeedOperator_Ref *impl) { 219 CeedInt ierr; 220 CeedInt dim, elemsize, size; 221 CeedElemRestriction Erestrict; 222 CeedEvalMode emode; 223 CeedBasis basis; 224 225 for (CeedInt i=0; i<numinputfields; i++) { 226 // Skip active input 227 if (skipactive) { 228 CeedVector vec; 229 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 230 if (vec == CEED_VECTOR_ACTIVE) 231 continue; 232 } 233 // Get elemsize, emode, size 234 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 235 CeedChk(ierr); 236 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 237 CeedChk(ierr); 238 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 239 CeedChk(ierr); 240 ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr); 241 // Basis action 242 switch(emode) { 243 case CEED_EVAL_NONE: 244 ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_HOST, 245 CEED_USE_POINTER, 246 &impl->edata[i][e*Q*size]); CeedChk(ierr); 247 break; 248 case CEED_EVAL_INTERP: 249 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 250 ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST, 251 CEED_USE_POINTER, 252 &impl->edata[i][e*elemsize*size]); 253 CeedChk(ierr); 254 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, 255 CEED_EVAL_INTERP, impl->evecsin[i], 256 impl->qvecsin[i]); CeedChk(ierr); 257 break; 258 case CEED_EVAL_GRAD: 259 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 260 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 261 ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST, 262 CEED_USE_POINTER, 263 &impl->edata[i][e*elemsize*size/dim]); 264 CeedChk(ierr); 265 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, 266 CEED_EVAL_GRAD, impl->evecsin[i], 267 impl->qvecsin[i]); CeedChk(ierr); 268 break; 269 case CEED_EVAL_WEIGHT: 270 break; // No action 271 // LCOV_EXCL_START 272 case CEED_EVAL_DIV: 273 case CEED_EVAL_CURL: { 274 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); 275 CeedChk(ierr); 276 Ceed ceed; 277 ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 278 return CeedError(ceed, 1, "Ceed evaluation mode not implemented"); 279 // LCOV_EXCL_STOP 280 } 281 } 282 } 283 return 0; 284 } 285 286 //------------------------------------------------------------------------------ 287 // Output Basis Action 288 //------------------------------------------------------------------------------ 289 static inline int CeedOperatorOutputBasis_Ref(CeedInt e, CeedInt Q, 290 CeedQFunctionField *qfoutputfields, CeedOperatorField *opoutputfields, 291 CeedInt numinputfields, CeedInt numoutputfields, CeedOperator op, 292 CeedOperator_Ref *impl) { 293 CeedInt ierr; 294 CeedInt dim, elemsize, size; 295 CeedElemRestriction Erestrict; 296 CeedEvalMode emode; 297 CeedBasis basis; 298 299 for (CeedInt i=0; i<numoutputfields; i++) { 300 // Get elemsize, emode, size 301 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 302 CeedChk(ierr); 303 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 304 CeedChk(ierr); 305 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 306 CeedChk(ierr); 307 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); CeedChk(ierr); 308 // Basis action 309 switch(emode) { 310 case CEED_EVAL_NONE: 311 break; // No action 312 case CEED_EVAL_INTERP: 313 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 314 CeedChk(ierr); 315 ierr = CeedVectorSetArray(impl->evecsout[i], CEED_MEM_HOST, 316 CEED_USE_POINTER, 317 &impl->edata[i + numinputfields][e*elemsize*size]); 318 CeedChk(ierr); 319 ierr = CeedBasisApply(basis, 1, CEED_TRANSPOSE, 320 CEED_EVAL_INTERP, impl->qvecsout[i], 321 impl->evecsout[i]); CeedChk(ierr); 322 break; 323 case CEED_EVAL_GRAD: 324 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 325 CeedChk(ierr); 326 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 327 ierr = CeedVectorSetArray(impl->evecsout[i], CEED_MEM_HOST, 328 CEED_USE_POINTER, 329 &impl->edata[i + numinputfields][e*elemsize*size/dim]); 330 CeedChk(ierr); 331 ierr = CeedBasisApply(basis, 1, CEED_TRANSPOSE, 332 CEED_EVAL_GRAD, impl->qvecsout[i], 333 impl->evecsout[i]); CeedChk(ierr); 334 break; 335 // LCOV_EXCL_START 336 case CEED_EVAL_WEIGHT: { 337 Ceed ceed; 338 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 339 return CeedError(ceed, 1, "CEED_EVAL_WEIGHT cannot be an output " 340 "evaluation mode"); 341 } 342 case CEED_EVAL_DIV: 343 case CEED_EVAL_CURL: { 344 Ceed ceed; 345 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 346 return CeedError(ceed, 1, "Ceed evaluation mode not implemented"); 347 // LCOV_EXCL_STOP 348 } 349 } 350 } 351 return 0; 352 } 353 354 //------------------------------------------------------------------------------ 355 // Restore Input Vectors 356 //------------------------------------------------------------------------------ 357 static inline int CeedOperatorRestoreInputs_Ref(CeedInt numinputfields, 358 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 359 const bool skipactive, CeedOperator_Ref *impl) { 360 CeedInt ierr; 361 CeedEvalMode emode; 362 363 for (CeedInt i=0; i<numinputfields; i++) { 364 // Skip active inputs 365 if (skipactive) { 366 CeedVector vec; 367 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 368 if (vec == CEED_VECTOR_ACTIVE) 369 continue; 370 } 371 // Restore input 372 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 373 CeedChk(ierr); 374 if (emode == CEED_EVAL_WEIGHT) { // Skip 375 } else { 376 ierr = CeedVectorRestoreArrayRead(impl->evecs[i], 377 (const CeedScalar **) &impl->edata[i]); 378 CeedChk(ierr); 379 } 380 } 381 return 0; 382 } 383 384 //------------------------------------------------------------------------------ 385 // Operator Apply 386 //------------------------------------------------------------------------------ 387 static int CeedOperatorApply_Ref(CeedOperator op, CeedVector invec, 388 CeedVector outvec, CeedRequest *request) { 389 int ierr; 390 CeedOperator_Ref *impl; 391 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 392 CeedQFunction qf; 393 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 394 CeedInt Q, numelements, numinputfields, numoutputfields, size; 395 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 396 ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr); 397 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 398 CeedChk(ierr); 399 CeedOperatorField *opinputfields, *opoutputfields; 400 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 401 CeedChk(ierr); 402 CeedQFunctionField *qfinputfields, *qfoutputfields; 403 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 404 CeedChk(ierr); 405 CeedEvalMode emode; 406 CeedVector vec; 407 CeedElemRestriction Erestrict; 408 409 // Setup 410 ierr = CeedOperatorSetup_Ref(op); CeedChk(ierr); 411 412 // Input Evecs and Restriction 413 ierr = CeedOperatorSetupInputs_Ref(numinputfields, qfinputfields, 414 opinputfields, invec, false, impl, 415 request); CeedChk(ierr); 416 417 // Output Evecs 418 for (CeedInt i=0; i<numoutputfields; i++) { 419 ierr = CeedVectorGetArray(impl->evecs[i+impl->numein], CEED_MEM_HOST, 420 &impl->edata[i + numinputfields]); CeedChk(ierr); 421 } 422 423 // Loop through elements 424 for (CeedInt e=0; e<numelements; e++) { 425 // Output pointers 426 for (CeedInt i=0; i<numoutputfields; i++) { 427 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 428 CeedChk(ierr); 429 if (emode == CEED_EVAL_NONE) { 430 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); 431 CeedChk(ierr); 432 ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_HOST, 433 CEED_USE_POINTER, 434 &impl->edata[i + numinputfields][e*Q*size]); 435 CeedChk(ierr); 436 } 437 } 438 439 // Input basis apply 440 ierr = CeedOperatorInputBasis_Ref(e, Q, qfinputfields, opinputfields, 441 numinputfields, false, impl); 442 CeedChk(ierr); 443 444 // Q function 445 if (!impl->identityqf) { 446 ierr = CeedQFunctionApply(qf, Q, impl->qvecsin, impl->qvecsout); 447 CeedChk(ierr); 448 } 449 450 // Output basis apply 451 ierr = CeedOperatorOutputBasis_Ref(e, Q, qfoutputfields, opoutputfields, 452 numinputfields, numoutputfields, op, impl); 453 CeedChk(ierr); 454 } 455 456 // Output restriction 457 for (CeedInt i=0; i<numoutputfields; i++) { 458 // Restore evec 459 ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein], 460 &impl->edata[i + numinputfields]); 461 CeedChk(ierr); 462 // Get output vector 463 ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr); 464 // Active 465 if (vec == CEED_VECTOR_ACTIVE) 466 vec = outvec; 467 // Restrict 468 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 469 CeedChk(ierr); 470 ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE, 471 impl->evecs[i+impl->numein], vec, request); 472 CeedChk(ierr); 473 } 474 475 // Restore input arrays 476 ierr = CeedOperatorRestoreInputs_Ref(numinputfields, qfinputfields, 477 opinputfields, false, impl); 478 CeedChk(ierr); 479 480 return 0; 481 } 482 483 //------------------------------------------------------------------------------ 484 // Assemble Linear QFunction 485 //------------------------------------------------------------------------------ 486 static int CeedOperatorLinearAssembleQFunction_Ref(CeedOperator op, 487 CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 488 int ierr; 489 CeedOperator_Ref *impl; 490 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 491 CeedQFunction qf; 492 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 493 CeedInt Q, numelements, numinputfields, numoutputfields, size; 494 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 495 ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr); 496 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 497 CeedChk(ierr); 498 CeedOperatorField *opinputfields, *opoutputfields; 499 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 500 CeedChk(ierr); 501 CeedQFunctionField *qfinputfields, *qfoutputfields; 502 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 503 CeedChk(ierr); 504 CeedVector vec; 505 CeedInt numactivein = 0, numactiveout = 0; 506 CeedVector *activein = NULL; 507 CeedScalar *a, *tmp; 508 Ceed ceed, ceedparent; 509 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 510 ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent); CeedChk(ierr); 511 ceedparent = ceedparent ? ceedparent : ceed; 512 513 // Setup 514 ierr = CeedOperatorSetup_Ref(op); CeedChk(ierr); 515 516 // Check for identity 517 if (impl->identityqf) 518 // LCOV_EXCL_START 519 return CeedError(ceed, 1, "Assembling identity QFunctions not supported"); 520 // LCOV_EXCL_STOP 521 522 // Input Evecs and Restriction 523 ierr = CeedOperatorSetupInputs_Ref(numinputfields, qfinputfields, 524 opinputfields, NULL, true, impl, request); 525 CeedChk(ierr); 526 527 // Count number of active input fields 528 for (CeedInt i=0; i<numinputfields; i++) { 529 // Get input vector 530 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 531 // Check if active input 532 if (vec == CEED_VECTOR_ACTIVE) { 533 ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr); 534 ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChk(ierr); 535 ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_HOST, &tmp); 536 CeedChk(ierr); 537 ierr = CeedRealloc(numactivein + size, &activein); CeedChk(ierr); 538 for (CeedInt field=0; field<size; field++) { 539 ierr = CeedVectorCreate(ceed, Q, &activein[numactivein+field]); 540 CeedChk(ierr); 541 ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_HOST, 542 CEED_USE_POINTER, &tmp[field*Q]); 543 CeedChk(ierr); 544 } 545 numactivein += size; 546 ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChk(ierr); 547 } 548 } 549 550 // Count number of active output fields 551 for (CeedInt i=0; i<numoutputfields; i++) { 552 // Get output vector 553 ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr); 554 // Check if active output 555 if (vec == CEED_VECTOR_ACTIVE) { 556 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); CeedChk(ierr); 557 numactiveout += size; 558 } 559 } 560 561 // Check sizes 562 if (!numactivein || !numactiveout) 563 // LCOV_EXCL_START 564 return CeedError(ceed, 1, "Cannot assemble QFunction without active inputs " 565 "and outputs"); 566 // LCOV_EXCL_STOP 567 568 // Create output restriction 569 CeedInt strides[3] = {1, Q, numactivein*numactiveout*Q}; /* *NOPAD* */ 570 ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q, 571 numactivein*numactiveout, 572 numactivein*numactiveout*numelements*Q, 573 strides, rstr); CeedChk(ierr); 574 // Create assembled vector 575 ierr = CeedVectorCreate(ceedparent, numelements*Q*numactivein*numactiveout, 576 assembled); CeedChk(ierr); 577 ierr = CeedVectorSetValue(*assembled, 0.0); CeedChk(ierr); 578 ierr = CeedVectorGetArray(*assembled, CEED_MEM_HOST, &a); CeedChk(ierr); 579 580 // Loop through elements 581 for (CeedInt e=0; e<numelements; e++) { 582 // Input basis apply 583 ierr = CeedOperatorInputBasis_Ref(e, Q, qfinputfields, opinputfields, 584 numinputfields, true, impl); 585 CeedChk(ierr); 586 587 // Assemble QFunction 588 for (CeedInt in=0; in<numactivein; in++) { 589 // Set Inputs 590 ierr = CeedVectorSetValue(activein[in], 1.0); CeedChk(ierr); 591 if (numactivein > 1) { 592 ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein], 593 0.0); CeedChk(ierr); 594 } 595 // Set Outputs 596 for (CeedInt out=0; out<numoutputfields; out++) { 597 // Get output vector 598 ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 599 CeedChk(ierr); 600 // Check if active output 601 if (vec == CEED_VECTOR_ACTIVE) { 602 CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST, 603 CEED_USE_POINTER, a); CeedChk(ierr); 604 ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size); 605 CeedChk(ierr); 606 a += size*Q; // Advance the pointer by the size of the output 607 } 608 } 609 // Apply QFunction 610 ierr = CeedQFunctionApply(qf, Q, impl->qvecsin, impl->qvecsout); 611 CeedChk(ierr); 612 } 613 } 614 615 // Un-set output Qvecs to prevent accidental overwrite of Assembled 616 for (CeedInt out=0; out<numoutputfields; out++) { 617 // Get output vector 618 ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 619 CeedChk(ierr); 620 // Check if active output 621 if (vec == CEED_VECTOR_ACTIVE) { 622 CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST, CEED_COPY_VALUES, 623 NULL); CeedChk(ierr); 624 } 625 } 626 627 // Restore input arrays 628 ierr = CeedOperatorRestoreInputs_Ref(numinputfields, qfinputfields, 629 opinputfields, true, impl); 630 CeedChk(ierr); 631 632 // Restore output 633 ierr = CeedVectorRestoreArray(*assembled, &a); CeedChk(ierr); 634 635 // Cleanup 636 for (CeedInt i=0; i<numactivein; i++) { 637 ierr = CeedVectorDestroy(&activein[i]); CeedChk(ierr); 638 } 639 ierr = CeedFree(&activein); CeedChk(ierr); 640 641 return 0; 642 } 643 644 //------------------------------------------------------------------------------ 645 // Get Basis Emode Pointer 646 //------------------------------------------------------------------------------ 647 static inline void CeedOperatorGetBasisPointer_Ref(const CeedScalar **basisptr, 648 CeedEvalMode emode, const CeedScalar *identity, const CeedScalar *interp, 649 const CeedScalar *grad) { 650 switch (emode) { 651 case CEED_EVAL_NONE: 652 *basisptr = identity; 653 break; 654 case CEED_EVAL_INTERP: 655 *basisptr = interp; 656 break; 657 case CEED_EVAL_GRAD: 658 *basisptr = grad; 659 break; 660 case CEED_EVAL_WEIGHT: 661 case CEED_EVAL_DIV: 662 case CEED_EVAL_CURL: 663 break; // Caught by QF Assembly 664 } 665 } 666 667 //------------------------------------------------------------------------------ 668 // Create point block restriction 669 //------------------------------------------------------------------------------ 670 static int CreatePBRestriction(CeedElemRestriction rstr, 671 CeedElemRestriction *pbRstr) { 672 int ierr; 673 Ceed ceed; 674 ierr = CeedElemRestrictionGetCeed(rstr, &ceed); CeedChk(ierr); 675 const CeedInt *offsets; 676 ierr = CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets); 677 CeedChk(ierr); 678 679 // Expand offsets 680 CeedInt nelem, ncomp, elemsize, compstride, max = 1, *pbOffsets; 681 ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChk(ierr); 682 ierr = CeedElemRestrictionGetNumComponents(rstr, &ncomp); CeedChk(ierr); 683 ierr = CeedElemRestrictionGetElementSize(rstr, &elemsize); CeedChk(ierr); 684 ierr = CeedElemRestrictionGetCompStride(rstr, &compstride); CeedChk(ierr); 685 CeedInt shift = ncomp; 686 if (compstride != 1) 687 shift *= ncomp; 688 ierr = CeedCalloc(nelem*elemsize, &pbOffsets); CeedChk(ierr); 689 for (CeedInt i = 0; i < nelem*elemsize; i++) { 690 pbOffsets[i] = offsets[i]*shift; 691 if (pbOffsets[i] > max) 692 max = pbOffsets[i]; 693 } 694 695 // Create new restriction 696 ierr = CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp*ncomp, 1, 697 max + ncomp*ncomp, CEED_MEM_HOST, 698 CEED_OWN_POINTER, pbOffsets, pbRstr); 699 CeedChk(ierr); 700 701 // Cleanup 702 ierr = CeedElemRestrictionRestoreOffsets(rstr, &offsets); CeedChk(ierr); 703 704 return 0; 705 } 706 707 //------------------------------------------------------------------------------ 708 // Assemble diagonal common code 709 //------------------------------------------------------------------------------ 710 static inline int CeedOperatorAssembleDiagonalCore_Ref(CeedOperator op, 711 CeedVector *assembled, CeedRequest *request, const bool pointBlock, 712 const bool createVector) { 713 int ierr; 714 Ceed ceed; 715 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 716 717 // Assemble QFunction 718 CeedQFunction qf; 719 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 720 CeedInt numinputfields, numoutputfields; 721 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 722 CeedChk(ierr); 723 CeedVector assembledqf; 724 CeedElemRestriction rstr; 725 ierr = CeedOperatorLinearAssembleQFunction(op, &assembledqf, &rstr, request); 726 CeedChk(ierr); 727 ierr = CeedElemRestrictionDestroy(&rstr); CeedChk(ierr); 728 CeedScalar maxnorm = 0; 729 ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm); CeedChk(ierr); 730 731 // Determine active input basis 732 CeedOperatorField *opfields; 733 CeedQFunctionField *qffields; 734 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 735 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 736 CeedInt numemodein = 0, ncomp, dim = 1; 737 CeedEvalMode *emodein = NULL; 738 CeedBasis basisin = NULL; 739 CeedElemRestriction rstrin = NULL; 740 for (CeedInt i=0; i<numinputfields; i++) { 741 CeedVector vec; 742 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 743 if (vec == CEED_VECTOR_ACTIVE) { 744 CeedElemRestriction rstr; 745 ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChk(ierr); 746 ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChk(ierr); 747 ierr = CeedBasisGetDimension(basisin, &dim); CeedChk(ierr); 748 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 749 CeedChk(ierr); 750 if (rstrin && rstrin != rstr) 751 // LCOV_EXCL_START 752 return CeedError(ceed, 1, 753 "Multi-field non-composite operator diagonal assembly not supported"); 754 // LCOV_EXCL_STOP 755 rstrin = rstr; 756 CeedEvalMode emode; 757 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); 758 CeedChk(ierr); 759 switch (emode) { 760 case CEED_EVAL_NONE: 761 case CEED_EVAL_INTERP: 762 ierr = CeedRealloc(numemodein + 1, &emodein); CeedChk(ierr); 763 emodein[numemodein] = emode; 764 numemodein += 1; 765 break; 766 case CEED_EVAL_GRAD: 767 ierr = CeedRealloc(numemodein + dim, &emodein); CeedChk(ierr); 768 for (CeedInt d=0; d<dim; d++) 769 emodein[numemodein+d] = emode; 770 numemodein += dim; 771 break; 772 case CEED_EVAL_WEIGHT: 773 case CEED_EVAL_DIV: 774 case CEED_EVAL_CURL: 775 break; // Caught by QF Assembly 776 } 777 } 778 } 779 780 // Determine active output basis 781 ierr = CeedOperatorGetFields(op, NULL, &opfields); CeedChk(ierr); 782 ierr = CeedQFunctionGetFields(qf, NULL, &qffields); CeedChk(ierr); 783 CeedInt numemodeout = 0; 784 CeedEvalMode *emodeout = NULL; 785 CeedBasis basisout = NULL; 786 CeedElemRestriction rstrout = NULL; 787 for (CeedInt i=0; i<numoutputfields; i++) { 788 CeedVector vec; 789 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 790 if (vec == CEED_VECTOR_ACTIVE) { 791 CeedElemRestriction rstr; 792 ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChk(ierr); 793 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 794 CeedChk(ierr); 795 if (rstrout && rstrout != rstr) 796 // LCOV_EXCL_START 797 return CeedError(ceed, 1, 798 "Multi-field non-composite operator diagonal assembly not supported"); 799 // LCOV_EXCL_STOP 800 rstrout = rstr; 801 CeedEvalMode emode; 802 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 803 switch (emode) { 804 case CEED_EVAL_NONE: 805 case CEED_EVAL_INTERP: 806 ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChk(ierr); 807 emodeout[numemodeout] = emode; 808 numemodeout += 1; 809 break; 810 case CEED_EVAL_GRAD: 811 ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChk(ierr); 812 for (CeedInt d=0; d<dim; d++) 813 emodeout[numemodeout+d] = emode; 814 numemodeout += dim; 815 break; 816 case CEED_EVAL_WEIGHT: 817 case CEED_EVAL_DIV: 818 case CEED_EVAL_CURL: 819 break; // Caught by QF Assembly 820 } 821 } 822 } 823 824 // Assemble point-block diagonal restriction, if needed 825 CeedElemRestriction diagrstr = rstrout; 826 if (pointBlock) { 827 ierr = CreatePBRestriction(rstrout, &diagrstr); CeedChk(ierr); 828 } 829 830 // Create diagonal vector 831 CeedVector elemdiag; 832 ierr = CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag); 833 CeedChk(ierr); 834 if (createVector) { 835 ierr = CeedElemRestrictionCreateVector(diagrstr, assembled, NULL); 836 CeedChk(ierr); 837 } 838 839 // Assemble element operator diagonals 840 CeedScalar *elemdiagarray, *assembledqfarray; 841 ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChk(ierr); 842 ierr = CeedVectorGetArray(elemdiag, CEED_MEM_HOST, &elemdiagarray); 843 CeedChk(ierr); 844 ierr = CeedVectorGetArray(assembledqf, CEED_MEM_HOST, &assembledqfarray); 845 CeedChk(ierr); 846 CeedInt nelem, nnodes, nqpts; 847 ierr = CeedElemRestrictionGetNumElements(diagrstr, &nelem); CeedChk(ierr); 848 ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChk(ierr); 849 ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChk(ierr); 850 // Basis matrices 851 const CeedScalar *interpin, *interpout, *gradin, *gradout; 852 CeedScalar *identity = NULL; 853 bool evalNone = false; 854 for (CeedInt i=0; i<numemodein; i++) 855 evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE); 856 for (CeedInt i=0; i<numemodeout; i++) 857 evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE); 858 if (evalNone) { 859 ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChk(ierr); 860 for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++) 861 identity[i*nnodes+i] = 1.0; 862 } 863 ierr = CeedBasisGetInterp(basisin, &interpin); CeedChk(ierr); 864 ierr = CeedBasisGetInterp(basisout, &interpout); CeedChk(ierr); 865 ierr = CeedBasisGetGrad(basisin, &gradin); CeedChk(ierr); 866 ierr = CeedBasisGetGrad(basisout, &gradout); CeedChk(ierr); 867 // Compute the diagonal of B^T D B 868 // Each element 869 for (CeedInt e=0; e<nelem; e++) { 870 CeedInt dout = -1; 871 // Each basis eval mode pair 872 for (CeedInt eout=0; eout<numemodeout; eout++) { 873 const CeedScalar *bt = NULL; 874 if (emodeout[eout] == CEED_EVAL_GRAD) 875 dout += 1; 876 CeedOperatorGetBasisPointer_Ref(&bt, emodeout[eout], identity, interpout, 877 &gradout[dout*nqpts*nnodes]); 878 CeedInt din = -1; 879 for (CeedInt ein=0; ein<numemodein; ein++) { 880 const CeedScalar *b = NULL; 881 if (emodein[ein] == CEED_EVAL_GRAD) 882 din += 1; 883 CeedOperatorGetBasisPointer_Ref(&b, emodein[ein], identity, interpin, 884 &gradin[din*nqpts*nnodes]); 885 // Each component 886 for (CeedInt compOut=0; compOut<ncomp; compOut++) 887 // Each qpoint/node pair 888 for (CeedInt q=0; q<nqpts; q++) 889 if (pointBlock) { 890 // Point Block Diagonal 891 for (CeedInt compIn=0; compIn<ncomp; compIn++) { 892 const CeedScalar qfvalue = 893 assembledqfarray[((((e*numemodein+ein)*ncomp+compIn)* 894 numemodeout+eout)*ncomp+compOut)*nqpts+q]; 895 if (fabs(qfvalue) > maxnorm*1e-12) 896 for (CeedInt n=0; n<nnodes; n++) 897 elemdiagarray[((e*ncomp+compOut)*ncomp+compIn)*nnodes+n] += 898 bt[q*nnodes+n] * qfvalue * b[q*nnodes+n]; 899 } 900 } else { 901 // Diagonal Only 902 const CeedScalar qfvalue = 903 assembledqfarray[((((e*numemodein+ein)*ncomp+compOut)* 904 numemodeout+eout)*ncomp+compOut)*nqpts+q]; 905 if (fabs(qfvalue) > maxnorm*1e-12) 906 for (CeedInt n=0; n<nnodes; n++) 907 elemdiagarray[(e*ncomp+compOut)*nnodes+n] += 908 bt[q*nnodes+n] * qfvalue * b[q*nnodes+n]; 909 } 910 } 911 } 912 } 913 ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChk(ierr); 914 ierr = CeedVectorRestoreArray(assembledqf, &assembledqfarray); CeedChk(ierr); 915 916 // Assemble local operator diagonal 917 if (createVector) { 918 ierr = CeedVectorSetValue(*assembled, 0.0); CeedChk(ierr); 919 } 920 ierr = CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag, 921 *assembled, request); CeedChk(ierr); 922 923 // Cleanup 924 if (pointBlock) { 925 ierr = CeedElemRestrictionDestroy(&diagrstr); CeedChk(ierr); 926 } 927 ierr = CeedVectorDestroy(&assembledqf); CeedChk(ierr); 928 ierr = CeedVectorDestroy(&elemdiag); CeedChk(ierr); 929 ierr = CeedFree(&emodein); CeedChk(ierr); 930 ierr = CeedFree(&emodeout); CeedChk(ierr); 931 ierr = CeedFree(&identity); CeedChk(ierr); 932 933 return 0; 934 } 935 936 //------------------------------------------------------------------------------ 937 // Assemble composite diagonal common code 938 //------------------------------------------------------------------------------ 939 static inline int CeedOperatorLinearAssembleDiagonalCompositeCore_Ref( 940 CeedOperator op, CeedVector *assembled, CeedRequest *request, 941 const bool pointBlock) { 942 int ierr; 943 CeedInt numSub; 944 CeedOperator *subOperators; 945 ierr = CeedOperatorGetNumSub(op, &numSub); CeedChk(ierr); 946 ierr = CeedOperatorGetSubList(op, &subOperators); CeedChk(ierr); 947 for (CeedInt i = 0; i < numSub; i++) { 948 ierr = CeedOperatorAssembleDiagonalCore_Ref(subOperators[i], assembled, 949 request, pointBlock, !i); CeedChk(ierr); 950 } 951 return 0; 952 } 953 954 //------------------------------------------------------------------------------ 955 // Assemble Linear Diagonal 956 //------------------------------------------------------------------------------ 957 static int CeedOperatorLinearAssembleDiagonal_Ref(CeedOperator op, 958 CeedVector *assembled, CeedRequest *request) { 959 int ierr; 960 bool isComposite; 961 ierr = CeedOperatorIsComposite(op, &isComposite); CeedChk(ierr); 962 if (isComposite) { 963 return CeedOperatorLinearAssembleDiagonalCompositeCore_Ref(op, assembled, 964 request, false); 965 } else { 966 return CeedOperatorAssembleDiagonalCore_Ref(op, assembled, request, false, 967 true); 968 } 969 } 970 971 //------------------------------------------------------------------------------ 972 // Assemble Linear Point Block Diagonal 973 //------------------------------------------------------------------------------ 974 static int CeedOperatorLinearAssemblePointBlockDiagonal_Ref(CeedOperator op, 975 CeedVector *assembled, CeedRequest *request) { 976 int ierr; 977 bool isComposite; 978 ierr = CeedOperatorIsComposite(op, &isComposite); CeedChk(ierr); 979 if (isComposite) { 980 return CeedOperatorLinearAssembleDiagonalCompositeCore_Ref(op, assembled, 981 request, true); 982 } else { 983 return CeedOperatorAssembleDiagonalCore_Ref(op, assembled, request, true, 984 true); 985 } 986 } 987 988 //------------------------------------------------------------------------------ 989 // Create FDM Element Inverse 990 //------------------------------------------------------------------------------ 991 int CeedOperatorCreateFDMElementInverse_Ref(CeedOperator op, 992 CeedOperator *fdminv, CeedRequest *request) { 993 int ierr; 994 Ceed ceed, ceedparent; 995 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 996 ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent); CeedChk(ierr); 997 ceedparent = ceedparent ? ceedparent : ceed; 998 CeedQFunction qf; 999 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 1000 1001 // Determine active input basis 1002 bool interp = false, grad = false; 1003 CeedBasis basis = NULL; 1004 CeedElemRestriction rstr = NULL; 1005 CeedOperatorField *opfields; 1006 CeedQFunctionField *qffields; 1007 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 1008 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 1009 CeedInt numinputfields; 1010 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, NULL); CeedChk(ierr); 1011 for (CeedInt i=0; i<numinputfields; i++) { 1012 CeedVector vec; 1013 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 1014 if (vec == CEED_VECTOR_ACTIVE) { 1015 CeedEvalMode emode; 1016 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 1017 interp = interp || emode == CEED_EVAL_INTERP; 1018 grad = grad || emode == CEED_EVAL_GRAD; 1019 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 1020 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 1021 CeedChk(ierr); 1022 } 1023 } 1024 if (!basis) 1025 // LCOV_EXCL_START 1026 return CeedError(ceed, 1, "No active field set"); 1027 // LCOV_EXCL_STOP 1028 CeedInt P1d, Q1d, elemsize, nqpts, dim, ncomp = 1, nelem = 1, lsize = 1; 1029 ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 1030 ierr = CeedBasisGetNumNodes(basis, &elemsize); CeedChk(ierr); 1031 ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 1032 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpts); CeedChk(ierr); 1033 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 1034 ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr); 1035 ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChk(ierr); 1036 ierr = CeedElemRestrictionGetLVectorSize(rstr, &lsize); CeedChk(ierr); 1037 1038 // Build and diagonalize 1D Mass and Laplacian 1039 bool tensorbasis; 1040 ierr = CeedBasisIsTensor(basis, &tensorbasis); CeedChk(ierr); 1041 if (!tensorbasis) 1042 // LCOV_EXCL_START 1043 return CeedError(ceed, 1, "FDMElementInverse only supported for tensor " 1044 "bases"); 1045 // LCOV_EXCL_STOP 1046 CeedScalar *work, *mass, *laplace, *x, *x2, *lambda; 1047 ierr = CeedMalloc(Q1d*P1d, &work); CeedChk(ierr); 1048 ierr = CeedMalloc(P1d*P1d, &mass); CeedChk(ierr); 1049 ierr = CeedMalloc(P1d*P1d, &laplace); CeedChk(ierr); 1050 ierr = CeedMalloc(P1d*P1d, &x); CeedChk(ierr); 1051 ierr = CeedMalloc(P1d*P1d, &x2); CeedChk(ierr); 1052 ierr = CeedMalloc(P1d, &lambda); CeedChk(ierr); 1053 // -- Mass 1054 const CeedScalar *interp1d, *grad1d, *qweight1d; 1055 ierr = CeedBasisGetInterp1D(basis, &interp1d); CeedChk(ierr); 1056 ierr = CeedBasisGetGrad1D(basis, &grad1d); CeedChk(ierr); 1057 ierr = CeedBasisGetQWeights(basis, &qweight1d); CeedChk(ierr); 1058 for (CeedInt i=0; i<Q1d; i++) 1059 for (CeedInt j=0; j<P1d; j++) 1060 work[i+j*Q1d] = interp1d[i*P1d+j]*qweight1d[i]; 1061 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)work, 1062 (const CeedScalar *)interp1d, mass, P1d, P1d, Q1d); 1063 CeedChk(ierr); 1064 // -- Laplacian 1065 for (CeedInt i=0; i<Q1d; i++) 1066 for (CeedInt j=0; j<P1d; j++) 1067 work[i+j*Q1d] = grad1d[i*P1d+j]*qweight1d[i]; 1068 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)work, 1069 (const CeedScalar *)grad1d, laplace, P1d, P1d, Q1d); 1070 CeedChk(ierr); 1071 // -- Diagonalize 1072 ierr = CeedSimultaneousDiagonalization(ceed, laplace, mass, x, lambda, P1d); 1073 CeedChk(ierr); 1074 ierr = CeedFree(&work); CeedChk(ierr); 1075 ierr = CeedFree(&mass); CeedChk(ierr); 1076 ierr = CeedFree(&laplace); CeedChk(ierr); 1077 for (CeedInt i=0; i<P1d; i++) 1078 for (CeedInt j=0; j<P1d; j++) 1079 x2[i+j*P1d] = x[j+i*P1d]; 1080 ierr = CeedFree(&x); CeedChk(ierr); 1081 1082 // Assemble QFunction 1083 CeedVector assembled; 1084 CeedElemRestriction rstr_qf; 1085 ierr = CeedOperatorLinearAssembleQFunction(op, &assembled, &rstr_qf, 1086 request); CeedChk(ierr); 1087 ierr = CeedElemRestrictionDestroy(&rstr_qf); CeedChk(ierr); 1088 CeedScalar maxnorm = 0; 1089 ierr = CeedVectorNorm(assembled, CEED_NORM_MAX, &maxnorm); CeedChk(ierr); 1090 1091 // Calculate element averages 1092 CeedInt nfields = ((interp?1:0) + (grad?dim:0))*((interp?1:0) + (grad?dim:0)); 1093 CeedScalar *elemavg; 1094 const CeedScalar *assembledarray, *qweightsarray; 1095 CeedVector qweights; 1096 ierr = CeedVectorCreate(ceedparent, nqpts, &qweights); CeedChk(ierr); 1097 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, 1098 CEED_VECTOR_NONE, qweights); CeedChk(ierr); 1099 ierr = CeedVectorGetArrayRead(assembled, CEED_MEM_HOST, &assembledarray); 1100 CeedChk(ierr); 1101 ierr = CeedVectorGetArrayRead(qweights, CEED_MEM_HOST, &qweightsarray); 1102 CeedChk(ierr); 1103 ierr = CeedCalloc(nelem, &elemavg); CeedChk(ierr); 1104 for (CeedInt e=0; e<nelem; e++) { 1105 CeedInt count = 0; 1106 for (CeedInt q=0; q<nqpts; q++) 1107 for (CeedInt i=0; i<ncomp*ncomp*nfields; i++) 1108 if (fabs(assembledarray[e*nelem*nqpts*ncomp*ncomp*nfields + 1109 i*nqpts + q]) > maxnorm*1e-12) { 1110 elemavg[e] += assembledarray[e*nelem*nqpts*ncomp*ncomp*nfields + 1111 i*nqpts + q] / qweightsarray[q]; 1112 count++; 1113 } 1114 if (count) 1115 elemavg[e] /= count; 1116 } 1117 ierr = CeedVectorRestoreArrayRead(assembled, &assembledarray); CeedChk(ierr); 1118 ierr = CeedVectorDestroy(&assembled); CeedChk(ierr); 1119 ierr = CeedVectorRestoreArrayRead(qweights, &qweightsarray); CeedChk(ierr); 1120 ierr = CeedVectorDestroy(&qweights); CeedChk(ierr); 1121 1122 // Build FDM diagonal 1123 CeedVector qdata; 1124 CeedScalar *qdataarray; 1125 ierr = CeedVectorCreate(ceedparent, nelem*ncomp*lsize, &qdata); CeedChk(ierr); 1126 ierr = CeedVectorSetArray(qdata, CEED_MEM_HOST, CEED_COPY_VALUES, NULL); 1127 CeedChk(ierr); 1128 ierr = CeedVectorGetArray(qdata, CEED_MEM_HOST, &qdataarray); CeedChk(ierr); 1129 for (CeedInt e=0; e<nelem; e++) 1130 for (CeedInt c=0; c<ncomp; c++) 1131 for (CeedInt n=0; n<lsize; n++) { 1132 if (interp) 1133 qdataarray[(e*ncomp+c)*lsize+n] = 1; 1134 if (grad) 1135 for (CeedInt d=0; d<dim; d++) { 1136 CeedInt i = (n / CeedIntPow(P1d, d)) % P1d; 1137 qdataarray[(e*ncomp+c)*lsize+n] += lambda[i]; 1138 } 1139 qdataarray[(e*ncomp+c)*lsize+n] = 1 / (elemavg[e] * 1140 qdataarray[(e*ncomp+c)*lsize+n]); 1141 } 1142 ierr = CeedFree(&elemavg); CeedChk(ierr); 1143 ierr = CeedVectorRestoreArray(qdata, &qdataarray); CeedChk(ierr); 1144 1145 // Setup FDM operator 1146 // -- Basis 1147 CeedBasis fdm_basis; 1148 CeedScalar *graddummy, *qrefdummy, *qweightdummy; 1149 ierr = CeedCalloc(P1d*P1d, &graddummy); CeedChk(ierr); 1150 ierr = CeedCalloc(P1d, &qrefdummy); CeedChk(ierr); 1151 ierr = CeedCalloc(P1d, &qweightdummy); CeedChk(ierr); 1152 ierr = CeedBasisCreateTensorH1(ceedparent, dim, ncomp, P1d, P1d, x2, 1153 graddummy, qrefdummy, qweightdummy, 1154 &fdm_basis); CeedChk(ierr); 1155 ierr = CeedFree(&graddummy); CeedChk(ierr); 1156 ierr = CeedFree(&qrefdummy); CeedChk(ierr); 1157 ierr = CeedFree(&qweightdummy); CeedChk(ierr); 1158 ierr = CeedFree(&x2); CeedChk(ierr); 1159 ierr = CeedFree(&lambda); CeedChk(ierr); 1160 1161 // -- Restriction 1162 CeedElemRestriction rstr_i; 1163 CeedInt strides[3] = {1, lsize, lsize*ncomp}; 1164 ierr = CeedElemRestrictionCreateStrided(ceedparent, nelem, lsize, ncomp, 1165 lsize*nelem*ncomp, strides, &rstr_i); 1166 CeedChk(ierr); 1167 // -- QFunction 1168 CeedQFunction mass_qf; 1169 ierr = CeedQFunctionCreateInteriorByName(ceedparent, "MassApply", &mass_qf); 1170 CeedChk(ierr); 1171 // -- Operator 1172 ierr = CeedOperatorCreate(ceedparent, mass_qf, NULL, NULL, fdminv); 1173 CeedChk(ierr); 1174 CeedOperatorSetField(*fdminv, "u", rstr_i, fdm_basis, CEED_VECTOR_ACTIVE); 1175 CeedChk(ierr); 1176 CeedOperatorSetField(*fdminv, "qdata", rstr_i, CEED_BASIS_COLLOCATED, qdata); 1177 CeedChk(ierr); 1178 CeedOperatorSetField(*fdminv, "v", rstr_i, fdm_basis, CEED_VECTOR_ACTIVE); 1179 CeedChk(ierr); 1180 1181 // Cleanup 1182 ierr = CeedVectorDestroy(&qdata); CeedChk(ierr); 1183 ierr = CeedBasisDestroy(&fdm_basis); CeedChk(ierr); 1184 ierr = CeedElemRestrictionDestroy(&rstr_i); CeedChk(ierr); 1185 ierr = CeedQFunctionDestroy(&mass_qf); CeedChk(ierr); 1186 1187 return 0; 1188 } 1189 1190 //------------------------------------------------------------------------------ 1191 // Operator Destroy 1192 //------------------------------------------------------------------------------ 1193 static int CeedOperatorDestroy_Ref(CeedOperator op) { 1194 int ierr; 1195 CeedOperator_Ref *impl; 1196 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 1197 1198 for (CeedInt i=0; i<impl->numein+impl->numeout; i++) { 1199 ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChk(ierr); 1200 } 1201 ierr = CeedFree(&impl->evecs); CeedChk(ierr); 1202 ierr = CeedFree(&impl->edata); CeedChk(ierr); 1203 ierr = CeedFree(&impl->inputstate); CeedChk(ierr); 1204 1205 for (CeedInt i=0; i<impl->numein; i++) { 1206 ierr = CeedVectorDestroy(&impl->evecsin[i]); CeedChk(ierr); 1207 ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChk(ierr); 1208 } 1209 ierr = CeedFree(&impl->evecsin); CeedChk(ierr); 1210 ierr = CeedFree(&impl->qvecsin); CeedChk(ierr); 1211 1212 for (CeedInt i=0; i<impl->numeout; i++) { 1213 ierr = CeedVectorDestroy(&impl->evecsout[i]); CeedChk(ierr); 1214 ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr); 1215 } 1216 ierr = CeedFree(&impl->evecsout); CeedChk(ierr); 1217 ierr = CeedFree(&impl->qvecsout); CeedChk(ierr); 1218 1219 ierr = CeedFree(&impl); CeedChk(ierr); 1220 return 0; 1221 } 1222 1223 //------------------------------------------------------------------------------ 1224 // Operator Create 1225 //------------------------------------------------------------------------------ 1226 int CeedOperatorCreate_Ref(CeedOperator op) { 1227 int ierr; 1228 Ceed ceed; 1229 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 1230 CeedOperator_Ref *impl; 1231 1232 ierr = CeedCalloc(1, &impl); CeedChk(ierr); 1233 ierr = CeedOperatorSetData(op, (void *)&impl); CeedChk(ierr); 1234 1235 ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", 1236 CeedOperatorLinearAssembleQFunction_Ref); 1237 CeedChk(ierr); 1238 ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleDiagonal", 1239 CeedOperatorLinearAssembleDiagonal_Ref); 1240 CeedChk(ierr); 1241 ierr = CeedSetBackendFunction(ceed, "Operator", op, 1242 "LinearAssemblePointBlockDiagonal", 1243 CeedOperatorLinearAssemblePointBlockDiagonal_Ref); 1244 CeedChk(ierr); 1245 ierr = CeedSetBackendFunction(ceed, "Operator", op, "CreateFDMElementInverse", 1246 CeedOperatorCreateFDMElementInverse_Ref); 1247 CeedChk(ierr); 1248 ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", 1249 CeedOperatorApply_Ref); CeedChk(ierr); 1250 ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy", 1251 CeedOperatorDestroy_Ref); CeedChk(ierr); 1252 return 0; 1253 } 1254 1255 //------------------------------------------------------------------------------ 1256 // Composite Operator Create 1257 //------------------------------------------------------------------------------ 1258 int CeedCompositeOperatorCreate_Ref(CeedOperator op) { 1259 int ierr; 1260 Ceed ceed; 1261 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 1262 ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleDiagonal", 1263 CeedOperatorLinearAssembleDiagonal_Ref); 1264 CeedChk(ierr); 1265 ierr = CeedSetBackendFunction(ceed, "Operator", op, 1266 "LinearAssemblePointBlockDiagonal", 1267 CeedOperatorLinearAssemblePointBlockDiagonal_Ref); 1268 CeedChk(ierr); 1269 return 0; 1270 } 1271 //------------------------------------------------------------------------------ 1272