1*5754ecacSJeremy L Thompson /// @file 2*5754ecacSJeremy L Thompson /// libCEED setup for solid mechanics example using PETSc 3*5754ecacSJeremy L Thompson 4*5754ecacSJeremy L Thompson #include "../include/setup-libceed.h" 5*5754ecacSJeremy L Thompson #include "../include/structs.h" 6*5754ecacSJeremy L Thompson #include "../include/utils.h" 7*5754ecacSJeremy L Thompson #include "../qfunctions/traction-boundary.h" // Traction boundaries 8*5754ecacSJeremy L Thompson #include "../qfunctions/constant-force.h" // Constant forcing function 9*5754ecacSJeremy L Thompson #include "../qfunctions/manufactured-force.h" // Manufactured solution forcing 10*5754ecacSJeremy L Thompson 11*5754ecacSJeremy L Thompson #if PETSC_VERSION_LT(3,14,0) 12*5754ecacSJeremy L Thompson # define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i) 13*5754ecacSJeremy L Thompson # define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i) 14*5754ecacSJeremy L Thompson #endif 15*5754ecacSJeremy L Thompson 16*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 17*5754ecacSJeremy L Thompson // Problem options 18*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 19*5754ecacSJeremy L Thompson // Forcing function data 20*5754ecacSJeremy L Thompson forcingData forcing_options[3] = { 21*5754ecacSJeremy L Thompson [FORCE_NONE] = { 22*5754ecacSJeremy L Thompson .setup_forcing = NULL, 23*5754ecacSJeremy L Thompson .setup_forcing_loc = NULL 24*5754ecacSJeremy L Thompson }, 25*5754ecacSJeremy L Thompson [FORCE_CONST] = { 26*5754ecacSJeremy L Thompson .setup_forcing = SetupConstantForce, 27*5754ecacSJeremy L Thompson .setup_forcing_loc = SetupConstantForce_loc 28*5754ecacSJeremy L Thompson }, 29*5754ecacSJeremy L Thompson [FORCE_MMS] = { 30*5754ecacSJeremy L Thompson .setup_forcing = SetupMMSForce, 31*5754ecacSJeremy L Thompson .setup_forcing_loc = SetupMMSForce_loc 32*5754ecacSJeremy L Thompson } 33*5754ecacSJeremy L Thompson }; 34*5754ecacSJeremy L Thompson 35*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 36*5754ecacSJeremy L Thompson // libCEED Functions 37*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 38*5754ecacSJeremy L Thompson // Destroy libCEED objects 39*5754ecacSJeremy L Thompson PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data) { 40*5754ecacSJeremy L Thompson PetscErrorCode ierr; 41*5754ecacSJeremy L Thompson 42*5754ecacSJeremy L Thompson PetscFunctionBegin; 43*5754ecacSJeremy L Thompson 44*5754ecacSJeremy L Thompson // Vectors 45*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->x_ceed); 46*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->y_ceed); 47*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->geo_data); 48*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < SOLIDS_MAX_NUMBER_FIELDS; i++) 49*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->stored_fields[i]); 50*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->geo_data_diagnostic); 51*5754ecacSJeremy L Thompson CeedVectorDestroy(&data->true_soln); 52*5754ecacSJeremy L Thompson // Restrictions 53*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_x); 54*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_u); 55*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_geo_data_i); 56*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < SOLIDS_MAX_NUMBER_FIELDS; i++) 57*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_stored_fields_i[i]); 58*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_energy); 59*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_diagnostic); 60*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&data->elem_restr_geo_data_diagnostic_i); 61*5754ecacSJeremy L Thompson // Bases 62*5754ecacSJeremy L Thompson CeedBasisDestroy(&data->basis_x); 63*5754ecacSJeremy L Thompson CeedBasisDestroy(&data->basis_u); 64*5754ecacSJeremy L Thompson CeedBasisDestroy(&data->basis_energy); 65*5754ecacSJeremy L Thompson CeedBasisDestroy(&data->basis_diagnostic); 66*5754ecacSJeremy L Thompson // QFunctions 67*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&data->qf_residual); 68*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&data->qf_jacobian); 69*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&data->qf_energy); 70*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&data->qf_diagnostic); 71*5754ecacSJeremy L Thompson // Operators 72*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_residual); 73*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_jacobian); 74*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_energy); 75*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_diagnostic); 76*5754ecacSJeremy L Thompson // Restriction and Prolongation data 77*5754ecacSJeremy L Thompson CeedBasisDestroy(&data->basis_c_to_f); 78*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_prolong); 79*5754ecacSJeremy L Thompson CeedOperatorDestroy(&data->op_restrict); 80*5754ecacSJeremy L Thompson 81*5754ecacSJeremy L Thompson ierr = PetscFree(data); CHKERRQ(ierr); 82*5754ecacSJeremy L Thompson 83*5754ecacSJeremy L Thompson PetscFunctionReturn(0); 84*5754ecacSJeremy L Thompson }; 85*5754ecacSJeremy L Thompson 86*5754ecacSJeremy L Thompson // Utility function - essential BC dofs are encoded in closure indices as -(i+1) 87*5754ecacSJeremy L Thompson PetscInt Involute(PetscInt i) { 88*5754ecacSJeremy L Thompson return i >= 0 ? i : -(i + 1); 89*5754ecacSJeremy L Thompson }; 90*5754ecacSJeremy L Thompson 91*5754ecacSJeremy L Thompson // Utility function to create local CEED restriction from DMPlex 92*5754ecacSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 93*5754ecacSJeremy L Thompson CeedInt height, DMLabel domain_label, CeedInt value, 94*5754ecacSJeremy L Thompson CeedElemRestriction *elem_restr) { 95*5754ecacSJeremy L Thompson 96*5754ecacSJeremy L Thompson PetscSection section; 97*5754ecacSJeremy L Thompson PetscInt p, num_elem, num_dof, *restr_indices, elem_offset, num_fields, dim, 98*5754ecacSJeremy L Thompson depth; 99*5754ecacSJeremy L Thompson DMLabel depth_label; 100*5754ecacSJeremy L Thompson IS depth_is, iter_is; 101*5754ecacSJeremy L Thompson Vec U_loc; 102*5754ecacSJeremy L Thompson const PetscInt *iter_indices; 103*5754ecacSJeremy L Thompson PetscErrorCode ierr; 104*5754ecacSJeremy L Thompson 105*5754ecacSJeremy L Thompson PetscFunctionBeginUser; 106*5754ecacSJeremy L Thompson 107*5754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 108*5754ecacSJeremy L Thompson dim -= height; 109*5754ecacSJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 110*5754ecacSJeremy L Thompson ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr); 111*5754ecacSJeremy L Thompson PetscInt num_comp[num_fields], field_offsets[num_fields+1]; 112*5754ecacSJeremy L Thompson field_offsets[0] = 0; 113*5754ecacSJeremy L Thompson for (PetscInt f = 0; f < num_fields; f++) { 114*5754ecacSJeremy L Thompson ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr); 115*5754ecacSJeremy L Thompson field_offsets[f+1] = field_offsets[f] + num_comp[f]; 116*5754ecacSJeremy L Thompson } 117*5754ecacSJeremy L Thompson 118*5754ecacSJeremy L Thompson ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 119*5754ecacSJeremy L Thompson ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr); 120*5754ecacSJeremy L Thompson ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is); 121*5754ecacSJeremy L Thompson CHKERRQ(ierr); 122*5754ecacSJeremy L Thompson if (domain_label) { 123*5754ecacSJeremy L Thompson IS domain_is; 124*5754ecacSJeremy L Thompson ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr); 125*5754ecacSJeremy L Thompson if (domain_is) { // domainIS is non-empty 126*5754ecacSJeremy L Thompson ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr); 127*5754ecacSJeremy L Thompson ierr = ISDestroy(&domain_is); CHKERRQ(ierr); 128*5754ecacSJeremy L Thompson } else { // domainIS is NULL (empty) 129*5754ecacSJeremy L Thompson iter_is = NULL; 130*5754ecacSJeremy L Thompson } 131*5754ecacSJeremy L Thompson ierr = ISDestroy(&depth_is); CHKERRQ(ierr); 132*5754ecacSJeremy L Thompson } else { 133*5754ecacSJeremy L Thompson iter_is = depth_is; 134*5754ecacSJeremy L Thompson } 135*5754ecacSJeremy L Thompson if (iter_is) { 136*5754ecacSJeremy L Thompson ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr); 137*5754ecacSJeremy L Thompson ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr); 138*5754ecacSJeremy L Thompson } else { 139*5754ecacSJeremy L Thompson num_elem = 0; 140*5754ecacSJeremy L Thompson iter_indices = NULL; 141*5754ecacSJeremy L Thompson } 142*5754ecacSJeremy L Thompson ierr = PetscMalloc1(num_elem*PetscPowInt(P, dim), &restr_indices); 143*5754ecacSJeremy L Thompson CHKERRQ(ierr); 144*5754ecacSJeremy L Thompson for (p = 0, elem_offset = 0; p < num_elem; p++) { 145*5754ecacSJeremy L Thompson PetscInt c = iter_indices[p]; 146*5754ecacSJeremy L Thompson PetscInt num_indices, *indices, num_nodes; 147*5754ecacSJeremy L Thompson ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 148*5754ecacSJeremy L Thompson &num_indices, &indices, NULL, NULL); 149*5754ecacSJeremy L Thompson CHKERRQ(ierr); 150*5754ecacSJeremy L Thompson bool flip = false; 151*5754ecacSJeremy L Thompson if (height > 0) { 152*5754ecacSJeremy L Thompson PetscInt num_cells, num_faces, start = -1; 153*5754ecacSJeremy L Thompson const PetscInt *orients, *faces, *cells; 154*5754ecacSJeremy L Thompson ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 155*5754ecacSJeremy L Thompson ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr); 156*5754ecacSJeremy L Thompson if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 157*5754ecacSJeremy L Thompson "Expected one cell in support of exterior face, but got %D cells", 158*5754ecacSJeremy L Thompson num_cells); 159*5754ecacSJeremy L Thompson ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 160*5754ecacSJeremy L Thompson ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr); 161*5754ecacSJeremy L Thompson for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;} 162*5754ecacSJeremy L Thompson if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, 163*5754ecacSJeremy L Thompson "Could not find face %D in cone of its support", 164*5754ecacSJeremy L Thompson c); 165*5754ecacSJeremy L Thompson ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 166*5754ecacSJeremy L Thompson if (orients[start] < 0) flip = true; 167*5754ecacSJeremy L Thompson } 168*5754ecacSJeremy L Thompson if (num_indices % field_offsets[num_fields]) SETERRQ1(PETSC_COMM_SELF, 169*5754ecacSJeremy L Thompson PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 170*5754ecacSJeremy L Thompson c); 171*5754ecacSJeremy L Thompson num_nodes = num_indices / field_offsets[num_fields]; 172*5754ecacSJeremy L Thompson for (PetscInt i = 0; i < num_nodes; i++) { 173*5754ecacSJeremy L Thompson PetscInt ii = i; 174*5754ecacSJeremy L Thompson if (flip) { 175*5754ecacSJeremy L Thompson if (P == num_nodes) ii = num_nodes - 1 - i; 176*5754ecacSJeremy L Thompson else if (P*P == num_nodes) { 177*5754ecacSJeremy L Thompson PetscInt row = i / P, col = i % P; 178*5754ecacSJeremy L Thompson ii = row + col * P; 179*5754ecacSJeremy L Thompson } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, 180*5754ecacSJeremy L Thompson "No support for flipping point with %D nodes != P (%D) or P^2", 181*5754ecacSJeremy L Thompson num_nodes, P); 182*5754ecacSJeremy L Thompson } 183*5754ecacSJeremy L Thompson // Check that indices are blocked by node and thus can be coalesced as a single field with 184*5754ecacSJeremy L Thompson // field_offsets[num_fields] = sum(num_comp) components. 185*5754ecacSJeremy L Thompson for (PetscInt f = 0; f < num_fields; f++) { 186*5754ecacSJeremy L Thompson for (PetscInt j = 0; j < num_comp[f]; j++) { 187*5754ecacSJeremy L Thompson if (Involute(indices[field_offsets[f]*num_nodes + ii*num_comp[f] + j]) 188*5754ecacSJeremy L Thompson != Involute(indices[ii*num_comp[0]]) + field_offsets[f] + j) 189*5754ecacSJeremy L Thompson SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 190*5754ecacSJeremy L Thompson "Cell %D closure indices not interlaced for node %D field %D component %D", 191*5754ecacSJeremy L Thompson c, ii, f, j); 192*5754ecacSJeremy L Thompson } 193*5754ecacSJeremy L Thompson } 194*5754ecacSJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 195*5754ecacSJeremy L Thompson PetscInt loc = Involute(indices[ii*num_comp[0]]); 196*5754ecacSJeremy L Thompson restr_indices[elem_offset++] = loc; 197*5754ecacSJeremy L Thompson } 198*5754ecacSJeremy L Thompson ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 199*5754ecacSJeremy L Thompson &num_indices, &indices, NULL, NULL); 200*5754ecacSJeremy L Thompson CHKERRQ(ierr); 201*5754ecacSJeremy L Thompson } 202*5754ecacSJeremy L Thompson if (elem_offset != num_elem*PetscPowInt(P, dim)) 203*5754ecacSJeremy L Thompson SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 204*5754ecacSJeremy L Thompson "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem, 205*5754ecacSJeremy L Thompson PetscPowInt(P, dim),elem_offset); 206*5754ecacSJeremy L Thompson if (iter_is) { 207*5754ecacSJeremy L Thompson ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr); 208*5754ecacSJeremy L Thompson } 209*5754ecacSJeremy L Thompson ierr = ISDestroy(&iter_is); CHKERRQ(ierr); 210*5754ecacSJeremy L Thompson 211*5754ecacSJeremy L Thompson ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr); 212*5754ecacSJeremy L Thompson ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr); 213*5754ecacSJeremy L Thompson ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr); 214*5754ecacSJeremy L Thompson CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, dim), 215*5754ecacSJeremy L Thompson field_offsets[num_fields], 216*5754ecacSJeremy L Thompson 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, restr_indices, 217*5754ecacSJeremy L Thompson elem_restr); 218*5754ecacSJeremy L Thompson ierr = PetscFree(restr_indices); CHKERRQ(ierr); 219*5754ecacSJeremy L Thompson PetscFunctionReturn(0); 220*5754ecacSJeremy L Thompson }; 221*5754ecacSJeremy L Thompson 222*5754ecacSJeremy L Thompson // Utility function to get Ceed Restriction for each domain 223*5754ecacSJeremy L Thompson PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, 224*5754ecacSJeremy L Thompson DMLabel domain_label, PetscInt value, 225*5754ecacSJeremy L Thompson CeedInt P, CeedInt Q, CeedInt q_data_size, 226*5754ecacSJeremy L Thompson CeedElemRestriction *elem_restr_q, 227*5754ecacSJeremy L Thompson CeedElemRestriction *elem_restr_x, 228*5754ecacSJeremy L Thompson CeedElemRestriction *elem_restr_qd_i) { 229*5754ecacSJeremy L Thompson 230*5754ecacSJeremy L Thompson DM dm_coord; 231*5754ecacSJeremy L Thompson CeedInt dim, num_local_elem; 232*5754ecacSJeremy L Thompson CeedInt Q_dim; 233*5754ecacSJeremy L Thompson PetscErrorCode ierr; 234*5754ecacSJeremy L Thompson 235*5754ecacSJeremy L Thompson PetscFunctionBeginUser; 236*5754ecacSJeremy L Thompson 237*5754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 238*5754ecacSJeremy L Thompson dim -= height; 239*5754ecacSJeremy L Thompson Q_dim = CeedIntPow(Q, dim); 240*5754ecacSJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 241*5754ecacSJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); 242*5754ecacSJeremy L Thompson CHKERRQ(ierr); 243*5754ecacSJeremy L Thompson if (elem_restr_q) { 244*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domain_label, value, 245*5754ecacSJeremy L Thompson elem_restr_q); CHKERRQ(ierr); 246*5754ecacSJeremy L Thompson } 247*5754ecacSJeremy L Thompson if (elem_restr_x) { 248*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm_coord, 2, height, domain_label, 249*5754ecacSJeremy L Thompson value, elem_restr_x); CHKERRQ(ierr); 250*5754ecacSJeremy L Thompson } 251*5754ecacSJeremy L Thompson if (elem_restr_qd_i) { 252*5754ecacSJeremy L Thompson CeedElemRestrictionGetNumElements(*elem_restr_q, &num_local_elem); 253*5754ecacSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_local_elem, Q_dim, 254*5754ecacSJeremy L Thompson q_data_size, q_data_size*num_local_elem*Q_dim, 255*5754ecacSJeremy L Thompson CEED_STRIDES_BACKEND, elem_restr_qd_i); 256*5754ecacSJeremy L Thompson } 257*5754ecacSJeremy L Thompson 258*5754ecacSJeremy L Thompson PetscFunctionReturn(0); 259*5754ecacSJeremy L Thompson }; 260*5754ecacSJeremy L Thompson 261*5754ecacSJeremy L Thompson // Set up libCEED on the fine grid for a given degree 262*5754ecacSJeremy L Thompson PetscErrorCode SetupLibceedFineLevel(DM dm, DM dm_energy, DM dm_diagnostic, 263*5754ecacSJeremy L Thompson Ceed ceed, AppCtx app_ctx, 264*5754ecacSJeremy L Thompson CeedQFunctionContext phys_ctx, 265*5754ecacSJeremy L Thompson ProblemData problem_data, 266*5754ecacSJeremy L Thompson PetscInt fine_level, PetscInt num_comp_u, 267*5754ecacSJeremy L Thompson PetscInt U_g_size, PetscInt U_loc_size, 268*5754ecacSJeremy L Thompson CeedVector force_ceed, 269*5754ecacSJeremy L Thompson CeedVector neumann_ceed, CeedData *data) { 270*5754ecacSJeremy L Thompson int ierr; 271*5754ecacSJeremy L Thompson CeedInt P = app_ctx->level_degrees[fine_level] + 1; 272*5754ecacSJeremy L Thompson CeedInt Q = app_ctx->level_degrees[fine_level] + 1 + app_ctx->q_extra; 273*5754ecacSJeremy L Thompson CeedInt dim, num_comp_x, num_comp_e = 1, num_comp_d = 5; 274*5754ecacSJeremy L Thompson CeedInt num_qpts; 275*5754ecacSJeremy L Thompson CeedInt geo_data_size = problem_data.geo_data_size; 276*5754ecacSJeremy L Thompson forcingType forcing_choice = app_ctx->forcing_choice; 277*5754ecacSJeremy L Thompson DM dm_coord; 278*5754ecacSJeremy L Thompson Vec coords; 279*5754ecacSJeremy L Thompson PetscInt c_start, c_end, num_elem; 280*5754ecacSJeremy L Thompson const PetscScalar *coordArray; 281*5754ecacSJeremy L Thompson CeedVector x_coord; 282*5754ecacSJeremy L Thompson CeedQFunction qf_setup_geo, qf_residual, qf_jacobian, qf_energy, qf_diagnostic; 283*5754ecacSJeremy L Thompson CeedOperator op_setup_geo, op_residual, op_jacobian, op_energy, op_diagnostic; 284*5754ecacSJeremy L Thompson 285*5754ecacSJeremy L Thompson PetscFunctionBeginUser; 286*5754ecacSJeremy L Thompson 287*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 288*5754ecacSJeremy L Thompson // libCEED bases 289*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 290*5754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 291*5754ecacSJeremy L Thompson num_comp_x = dim; 292*5754ecacSJeremy L Thompson // -- Coordinate basis 293*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, 294*5754ecacSJeremy L Thompson problem_data.quadrature_mode, 295*5754ecacSJeremy L Thompson &data[fine_level]->basis_x); 296*5754ecacSJeremy L Thompson // -- Solution basis 297*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, 298*5754ecacSJeremy L Thompson problem_data.quadrature_mode, 299*5754ecacSJeremy L Thompson &data[fine_level]->basis_u); 300*5754ecacSJeremy L Thompson // -- Energy basis 301*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_e, P, Q, 302*5754ecacSJeremy L Thompson problem_data.quadrature_mode, 303*5754ecacSJeremy L Thompson &data[fine_level]->basis_energy); 304*5754ecacSJeremy L Thompson // -- Diagnostic output basis 305*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, P, CEED_GAUSS_LOBATTO, 306*5754ecacSJeremy L Thompson &data[fine_level]->basis_diagnostic); 307*5754ecacSJeremy L Thompson 308*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 309*5754ecacSJeremy L Thompson // libCEED restrictions 310*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 311*5754ecacSJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 312*5754ecacSJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); 313*5754ecacSJeremy L Thompson CHKERRQ(ierr); 314*5754ecacSJeremy L Thompson 315*5754ecacSJeremy L Thompson // -- Coordinate restriction 316*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm_coord, 2, 0, 0, 0, 317*5754ecacSJeremy L Thompson &(data[fine_level]->elem_restr_x)); 318*5754ecacSJeremy L Thompson CHKERRQ(ierr); 319*5754ecacSJeremy L Thompson // -- Solution restriction 320*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, P, 0, 0, 0, 321*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_u); 322*5754ecacSJeremy L Thompson CHKERRQ(ierr); 323*5754ecacSJeremy L Thompson // -- Energy restriction 324*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm_energy, P, 0, 0, 0, 325*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_energy); 326*5754ecacSJeremy L Thompson CHKERRQ(ierr); 327*5754ecacSJeremy L Thompson // -- Diagnostic data restriction 328*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm_diagnostic, P, 0, 0, 0, 329*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_diagnostic); 330*5754ecacSJeremy L Thompson CHKERRQ(ierr); 331*5754ecacSJeremy L Thompson 332*5754ecacSJeremy L Thompson // -- Stored data at quadrature points 333*5754ecacSJeremy L Thompson ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr); 334*5754ecacSJeremy L Thompson num_elem = c_end - c_start; 335*5754ecacSJeremy L Thompson CeedBasisGetNumQuadraturePoints(data[fine_level]->basis_u, &num_qpts); 336*5754ecacSJeremy L Thompson // ---- Geometric data restriction, residual and Jacobian operators 337*5754ecacSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, geo_data_size, 338*5754ecacSJeremy L Thompson num_elem*num_qpts*geo_data_size, 339*5754ecacSJeremy L Thompson CEED_STRIDES_BACKEND, 340*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_geo_data_i); 341*5754ecacSJeremy L Thompson // ---- Stored field restrictions 342*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 343*5754ecacSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, 344*5754ecacSJeremy L Thompson problem_data.field_sizes[i], 345*5754ecacSJeremy L Thompson num_elem*num_qpts*problem_data.field_sizes[i], 346*5754ecacSJeremy L Thompson CEED_STRIDES_BACKEND, 347*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_stored_fields_i[i]); 348*5754ecacSJeremy L Thompson } 349*5754ecacSJeremy L Thompson // ---- Geometric data restriction, diagnostic operator 350*5754ecacSJeremy L Thompson CeedElemRestrictionCreateStrided(ceed, num_elem, P*P*P, geo_data_size, 351*5754ecacSJeremy L Thompson num_elem*P*P*P*geo_data_size, 352*5754ecacSJeremy L Thompson CEED_STRIDES_BACKEND, 353*5754ecacSJeremy L Thompson &data[fine_level]->elem_restr_geo_data_diagnostic_i); 354*5754ecacSJeremy L Thompson 355*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 356*5754ecacSJeremy L Thompson // Element coordinates 357*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 358*5754ecacSJeremy L Thompson ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr); 359*5754ecacSJeremy L Thompson ierr = VecGetArrayRead(coords, &coordArray); CHKERRQ(ierr); 360*5754ecacSJeremy L Thompson 361*5754ecacSJeremy L Thompson CeedElemRestrictionCreateVector(data[fine_level]->elem_restr_x, &x_coord, NULL); 362*5754ecacSJeremy L Thompson CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, 363*5754ecacSJeremy L Thompson (PetscScalar *)coordArray); 364*5754ecacSJeremy L Thompson ierr = VecRestoreArrayRead(coords, &coordArray); CHKERRQ(ierr); 365*5754ecacSJeremy L Thompson 366*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 367*5754ecacSJeremy L Thompson // Persistent libCEED vectors 368*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 369*5754ecacSJeremy L Thompson // -- Operator action variables 370*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, U_loc_size, &data[fine_level]->x_ceed); 371*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, U_loc_size, &data[fine_level]->y_ceed); 372*5754ecacSJeremy L Thompson // -- Geometric data vector 373*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, num_elem*num_qpts*geo_data_size, 374*5754ecacSJeremy L Thompson &data[fine_level]->geo_data); 375*5754ecacSJeremy L Thompson // -- Stored field vectors 376*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 377*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, num_elem*num_qpts*problem_data.field_sizes[i], 378*5754ecacSJeremy L Thompson &data[fine_level]->stored_fields[i]); 379*5754ecacSJeremy L Thompson } 380*5754ecacSJeremy L Thompson // -- Collocated geometric data vector 381*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, num_elem*P*P*P*geo_data_size, 382*5754ecacSJeremy L Thompson &data[fine_level]->geo_data_diagnostic); 383*5754ecacSJeremy L Thompson 384*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 385*5754ecacSJeremy L Thompson // Geometric factor computation 386*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 387*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the quadrature data 388*5754ecacSJeremy L Thompson // geo_data returns dXdx_i,j and w * det. 389*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 390*5754ecacSJeremy L Thompson // -- QFunction 391*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.setup_geo, 392*5754ecacSJeremy L Thompson problem_data.setup_geo_loc, &qf_setup_geo); 393*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*dim, CEED_EVAL_GRAD); 394*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT); 395*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_setup_geo, "geo_data", geo_data_size, CEED_EVAL_NONE); 396*5754ecacSJeremy L Thompson // -- Operator 397*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, 398*5754ecacSJeremy L Thompson CEED_QFUNCTION_NONE, &op_setup_geo); 399*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "dx", data[fine_level]->elem_restr_x, 400*5754ecacSJeremy L Thompson data[fine_level]->basis_x, CEED_VECTOR_ACTIVE); 401*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, 402*5754ecacSJeremy L Thompson data[fine_level]->basis_x, CEED_VECTOR_NONE); 403*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "geo_data", 404*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_i, 405*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 406*5754ecacSJeremy L Thompson // -- Compute the quadrature data 407*5754ecacSJeremy L Thompson CeedOperatorApply(op_setup_geo, x_coord, data[fine_level]->geo_data, 408*5754ecacSJeremy L Thompson CEED_REQUEST_IMMEDIATE); 409*5754ecacSJeremy L Thompson // -- Cleanup 410*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&qf_setup_geo); 411*5754ecacSJeremy L Thompson CeedOperatorDestroy(&op_setup_geo); 412*5754ecacSJeremy L Thompson 413*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 414*5754ecacSJeremy L Thompson // Local residual evaluator 415*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 416*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the residual of the 417*5754ecacSJeremy L Thompson // non-linear PDE. 418*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 419*5754ecacSJeremy L Thompson // -- QFunction 420*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.residual, 421*5754ecacSJeremy L Thompson problem_data.residual_loc, &qf_residual); 422*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_residual, "du", num_comp_u*dim, CEED_EVAL_GRAD); 423*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_residual, "geo_data", geo_data_size, CEED_EVAL_NONE); 424*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_residual, "dv", num_comp_u*dim, CEED_EVAL_GRAD); 425*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 426*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_residual, problem_data.field_names[i], 427*5754ecacSJeremy L Thompson problem_data.field_sizes[i], CEED_EVAL_NONE); 428*5754ecacSJeremy L Thompson } 429*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_residual, phys_ctx); 430*5754ecacSJeremy L Thompson // -- Operator 431*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_residual, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 432*5754ecacSJeremy L Thompson &op_residual); 433*5754ecacSJeremy L Thompson CeedOperatorSetField(op_residual, "du", data[fine_level]->elem_restr_u, 434*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 435*5754ecacSJeremy L Thompson CeedOperatorSetField(op_residual, "geo_data", 436*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_i, 437*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, data[fine_level]->geo_data); 438*5754ecacSJeremy L Thompson CeedOperatorSetField(op_residual, "dv", data[fine_level]->elem_restr_u, 439*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 440*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 441*5754ecacSJeremy L Thompson CeedOperatorSetField(op_residual, problem_data.field_names[i], 442*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_stored_fields_i[i], 443*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, 444*5754ecacSJeremy L Thompson data[fine_level]->stored_fields[i]); 445*5754ecacSJeremy L Thompson } 446*5754ecacSJeremy L Thompson // -- Save libCEED data 447*5754ecacSJeremy L Thompson data[fine_level]->qf_residual = qf_residual; 448*5754ecacSJeremy L Thompson data[fine_level]->op_residual = op_residual; 449*5754ecacSJeremy L Thompson 450*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 451*5754ecacSJeremy L Thompson // Jacobian evaluator 452*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 453*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the action of the 454*5754ecacSJeremy L Thompson // Jacobian for each linear solve. 455*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 456*5754ecacSJeremy L Thompson // -- QFunction 457*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.jacobian, 458*5754ecacSJeremy L Thompson problem_data.jacobian_loc, &qf_jacobian); 459*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_jacobian, "deltadu", num_comp_u*dim, CEED_EVAL_GRAD); 460*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_jacobian, "geo_data", geo_data_size, CEED_EVAL_NONE); 461*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 462*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_jacobian, problem_data.field_names[i], 463*5754ecacSJeremy L Thompson problem_data.field_sizes[i], CEED_EVAL_NONE); 464*5754ecacSJeremy L Thompson } 465*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_jacobian, "deltadv", num_comp_u*dim, CEED_EVAL_GRAD); 466*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_jacobian, phys_ctx); 467*5754ecacSJeremy L Thompson // -- Operator 468*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_jacobian, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 469*5754ecacSJeremy L Thompson &op_jacobian); 470*5754ecacSJeremy L Thompson CeedOperatorSetField(op_jacobian, "deltadu", data[fine_level]->elem_restr_u, 471*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 472*5754ecacSJeremy L Thompson CeedOperatorSetField(op_jacobian, "geo_data", 473*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_i, 474*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, data[fine_level]->geo_data); 475*5754ecacSJeremy L Thompson CeedOperatorSetField(op_jacobian, "deltadv", data[fine_level]->elem_restr_u, 476*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 477*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < problem_data.number_fields_stored; i++) { 478*5754ecacSJeremy L Thompson CeedOperatorSetField(op_jacobian, problem_data.field_names[i], 479*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_stored_fields_i[i], 480*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, 481*5754ecacSJeremy L Thompson data[fine_level]->stored_fields[i]); 482*5754ecacSJeremy L Thompson } 483*5754ecacSJeremy L Thompson // -- Save libCEED data 484*5754ecacSJeremy L Thompson data[fine_level]->qf_jacobian = qf_jacobian; 485*5754ecacSJeremy L Thompson data[fine_level]->op_jacobian = op_jacobian; 486*5754ecacSJeremy L Thompson 487*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 488*5754ecacSJeremy L Thompson // Traction boundary conditions, if needed 489*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 490*5754ecacSJeremy L Thompson if (app_ctx->bc_traction_count > 0) { 491*5754ecacSJeremy L Thompson // -- Setup 492*5754ecacSJeremy L Thompson DMLabel domain_label; 493*5754ecacSJeremy L Thompson ierr = DMGetLabel(dm, "Face Sets", &domain_label); CHKERRQ(ierr); 494*5754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 495*5754ecacSJeremy L Thompson 496*5754ecacSJeremy L Thompson // -- Basis 497*5754ecacSJeremy L Thompson CeedInt height = 1; 498*5754ecacSJeremy L Thompson CeedBasis basis_x_face, basis_u_face; 499*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim - height, num_comp_x, 2, Q, 500*5754ecacSJeremy L Thompson problem_data.quadrature_mode, &basis_x_face); 501*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim - height, num_comp_u, P, Q, 502*5754ecacSJeremy L Thompson problem_data.quadrature_mode, &basis_u_face); 503*5754ecacSJeremy L Thompson // -- QFunction 504*5754ecacSJeremy L Thompson CeedQFunction qf_traction; 505*5754ecacSJeremy L Thompson CeedQFunctionContext traction_ctx; 506*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, SetupTractionBCs, SetupTractionBCs_loc, 507*5754ecacSJeremy L Thompson &qf_traction); 508*5754ecacSJeremy L Thompson CeedQFunctionContextCreate(ceed, &traction_ctx); 509*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_traction, traction_ctx); 510*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_traction, "dx", num_comp_x*(num_comp_x - height), 511*5754ecacSJeremy L Thompson CEED_EVAL_GRAD); 512*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_traction, "weight", 1, CEED_EVAL_WEIGHT); 513*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_traction, "v", num_comp_u, CEED_EVAL_INTERP); 514*5754ecacSJeremy L Thompson 515*5754ecacSJeremy L Thompson // -- Compute contribution on each boundary face 516*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < app_ctx->bc_traction_count; i++) { 517*5754ecacSJeremy L Thompson CeedElemRestriction elem_restr_x_face, elem_restr_u_face; 518*5754ecacSJeremy L Thompson CeedOperator op_traction; 519*5754ecacSJeremy L Thompson CeedQFunctionContextSetData(traction_ctx, CEED_MEM_HOST, CEED_USE_POINTER, 520*5754ecacSJeremy L Thompson 3 * sizeof(CeedScalar), 521*5754ecacSJeremy L Thompson app_ctx->bc_traction_vector[i]); 522*5754ecacSJeremy L Thompson // Setup restriction 523*5754ecacSJeremy L Thompson ierr = GetRestrictionForDomain(ceed, dm, height, domain_label, 524*5754ecacSJeremy L Thompson app_ctx->bc_traction_faces[i], P, Q, 525*5754ecacSJeremy L Thompson 0, &elem_restr_u_face, &elem_restr_x_face, NULL); 526*5754ecacSJeremy L Thompson CHKERRQ(ierr); 527*5754ecacSJeremy L Thompson // ---- Create boundary Operator 528*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_traction, NULL, NULL, &op_traction); 529*5754ecacSJeremy L Thompson CeedOperatorSetField(op_traction, "dx", elem_restr_x_face, basis_x_face, 530*5754ecacSJeremy L Thompson CEED_VECTOR_ACTIVE); 531*5754ecacSJeremy L Thompson CeedOperatorSetField(op_traction, "weight", CEED_ELEMRESTRICTION_NONE, 532*5754ecacSJeremy L Thompson basis_x_face, CEED_VECTOR_NONE); 533*5754ecacSJeremy L Thompson CeedOperatorSetField(op_traction, "v", elem_restr_u_face, 534*5754ecacSJeremy L Thompson basis_u_face, CEED_VECTOR_ACTIVE); 535*5754ecacSJeremy L Thompson // ---- Compute traction on face 536*5754ecacSJeremy L Thompson CeedOperatorApplyAdd(op_traction, x_coord, neumann_ceed, 537*5754ecacSJeremy L Thompson CEED_REQUEST_IMMEDIATE); 538*5754ecacSJeremy L Thompson // ---- Cleanup 539*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&elem_restr_x_face); 540*5754ecacSJeremy L Thompson CeedElemRestrictionDestroy(&elem_restr_u_face); 541*5754ecacSJeremy L Thompson CeedOperatorDestroy(&op_traction); 542*5754ecacSJeremy L Thompson } 543*5754ecacSJeremy L Thompson // -- Cleanup 544*5754ecacSJeremy L Thompson CeedBasisDestroy(&basis_x_face); 545*5754ecacSJeremy L Thompson CeedBasisDestroy(&basis_u_face); 546*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&qf_traction); 547*5754ecacSJeremy L Thompson CeedQFunctionContextDestroy(&traction_ctx); 548*5754ecacSJeremy L Thompson } 549*5754ecacSJeremy L Thompson 550*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 551*5754ecacSJeremy L Thompson // Forcing term, if needed 552*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 553*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the forcing term (RHS) 554*5754ecacSJeremy L Thompson // for the non-linear PDE. 555*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 556*5754ecacSJeremy L Thompson if (forcing_choice != FORCE_NONE) { 557*5754ecacSJeremy L Thompson CeedQFunction qf_setup_force; 558*5754ecacSJeremy L Thompson CeedOperator op_setup_force; 559*5754ecacSJeremy L Thompson 560*5754ecacSJeremy L Thompson // -- QFunction 561*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, 562*5754ecacSJeremy L Thompson forcing_options[forcing_choice].setup_forcing, 563*5754ecacSJeremy L Thompson forcing_options[forcing_choice].setup_forcing_loc, 564*5754ecacSJeremy L Thompson &qf_setup_force); 565*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_force, "x", num_comp_x, CEED_EVAL_INTERP); 566*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_force, "geo_data", geo_data_size, 567*5754ecacSJeremy L Thompson CEED_EVAL_NONE); 568*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_setup_force, "force", num_comp_u, CEED_EVAL_INTERP); 569*5754ecacSJeremy L Thompson if (forcing_choice == FORCE_MMS) { 570*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_setup_force, phys_ctx); 571*5754ecacSJeremy L Thompson } else { 572*5754ecacSJeremy L Thompson CeedQFunctionContext ctxForcing; 573*5754ecacSJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxForcing); 574*5754ecacSJeremy L Thompson CeedQFunctionContextSetData(ctxForcing, CEED_MEM_HOST, CEED_USE_POINTER, 575*5754ecacSJeremy L Thompson sizeof(*app_ctx->forcing_vector), 576*5754ecacSJeremy L Thompson app_ctx->forcing_vector); 577*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_setup_force, ctxForcing); 578*5754ecacSJeremy L Thompson CeedQFunctionContextDestroy(&ctxForcing); 579*5754ecacSJeremy L Thompson } 580*5754ecacSJeremy L Thompson // -- Operator 581*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_setup_force, CEED_QFUNCTION_NONE, 582*5754ecacSJeremy L Thompson CEED_QFUNCTION_NONE, &op_setup_force); 583*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_force, "x", data[fine_level]->elem_restr_x, 584*5754ecacSJeremy L Thompson data[fine_level]->basis_x, CEED_VECTOR_ACTIVE); 585*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_force, "geo_data", 586*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_i, 587*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, data[fine_level]->geo_data); 588*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_force, "force", data[fine_level]->elem_restr_u, 589*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 590*5754ecacSJeremy L Thompson // -- Compute forcing term 591*5754ecacSJeremy L Thompson CeedOperatorApply(op_setup_force, x_coord, force_ceed, CEED_REQUEST_IMMEDIATE); 592*5754ecacSJeremy L Thompson // -- Cleanup 593*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&qf_setup_force); 594*5754ecacSJeremy L Thompson CeedOperatorDestroy(&op_setup_force); 595*5754ecacSJeremy L Thompson } 596*5754ecacSJeremy L Thompson 597*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 598*5754ecacSJeremy L Thompson // True solution, for MMS 599*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 600*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the true solution at 601*5754ecacSJeremy L Thompson // the mesh nodes for validation with the manufactured solution. 602*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 603*5754ecacSJeremy L Thompson if (problem_data.true_soln) { 604*5754ecacSJeremy L Thompson CeedScalar *true_array; 605*5754ecacSJeremy L Thompson const CeedScalar *mult_array; 606*5754ecacSJeremy L Thompson CeedVector mult_vec; 607*5754ecacSJeremy L Thompson CeedBasis basis_x_true; 608*5754ecacSJeremy L Thompson CeedQFunction qf_true; 609*5754ecacSJeremy L Thompson CeedOperator op_true; 610*5754ecacSJeremy L Thompson 611*5754ecacSJeremy L Thompson // -- Solution vector 612*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, U_loc_size, &(data[fine_level]->true_soln)); 613*5754ecacSJeremy L Thompson // -- Basis 614*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, CEED_GAUSS_LOBATTO, 615*5754ecacSJeremy L Thompson &basis_x_true); 616*5754ecacSJeremy L Thompson // QFunction 617*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.true_soln, 618*5754ecacSJeremy L Thompson problem_data.true_soln_loc, &qf_true); 619*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_true, "x", num_comp_x, CEED_EVAL_INTERP); 620*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_true, "true_soln", num_comp_u, CEED_EVAL_NONE); 621*5754ecacSJeremy L Thompson // Operator 622*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_true, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 623*5754ecacSJeremy L Thompson &op_true); 624*5754ecacSJeremy L Thompson CeedOperatorSetField(op_true, "x", data[fine_level]->elem_restr_x, basis_x_true, 625*5754ecacSJeremy L Thompson CEED_VECTOR_ACTIVE); 626*5754ecacSJeremy L Thompson CeedOperatorSetField(op_true, "true_soln", data[fine_level]->elem_restr_u, 627*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 628*5754ecacSJeremy L Thompson // -- Compute true solution 629*5754ecacSJeremy L Thompson CeedOperatorApply(op_true, x_coord, data[fine_level]->true_soln, 630*5754ecacSJeremy L Thompson CEED_REQUEST_IMMEDIATE); 631*5754ecacSJeremy L Thompson // -- Multiplicity calculation 632*5754ecacSJeremy L Thompson CeedElemRestrictionCreateVector(data[fine_level]->elem_restr_u, &mult_vec, 633*5754ecacSJeremy L Thompson NULL); 634*5754ecacSJeremy L Thompson CeedVectorSetValue(mult_vec, 0.); 635*5754ecacSJeremy L Thompson CeedElemRestrictionGetMultiplicity(data[fine_level]->elem_restr_u, mult_vec); 636*5754ecacSJeremy L Thompson // -- Multiplicity correction 637*5754ecacSJeremy L Thompson CeedVectorGetArray(data[fine_level]->true_soln, CEED_MEM_HOST, &true_array); 638*5754ecacSJeremy L Thompson CeedVectorGetArrayRead(mult_vec, CEED_MEM_HOST, &mult_array); 639*5754ecacSJeremy L Thompson for (CeedInt i = 0; i < U_loc_size; i++) 640*5754ecacSJeremy L Thompson true_array[i] /= mult_array[i]; 641*5754ecacSJeremy L Thompson CeedVectorRestoreArray(data[fine_level]->true_soln, &true_array); 642*5754ecacSJeremy L Thompson CeedVectorRestoreArrayRead(mult_vec, &mult_array); 643*5754ecacSJeremy L Thompson // -- Cleanup 644*5754ecacSJeremy L Thompson CeedVectorDestroy(&mult_vec); 645*5754ecacSJeremy L Thompson CeedBasisDestroy(&basis_x_true); 646*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&qf_true); 647*5754ecacSJeremy L Thompson CeedOperatorDestroy(&op_true); 648*5754ecacSJeremy L Thompson } 649*5754ecacSJeremy L Thompson 650*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 651*5754ecacSJeremy L Thompson // Local energy computation 652*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 653*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes the strain energy 654*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 655*5754ecacSJeremy L Thompson // -- QFunction 656*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.energy, 657*5754ecacSJeremy L Thompson problem_data.energy_loc, &qf_energy); 658*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_energy, "du", num_comp_u*dim, CEED_EVAL_GRAD); 659*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_energy, "geo_data", geo_data_size, CEED_EVAL_NONE); 660*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_energy, "energy", num_comp_e, CEED_EVAL_INTERP); 661*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_energy, phys_ctx); 662*5754ecacSJeremy L Thompson // -- Operator 663*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_energy, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 664*5754ecacSJeremy L Thompson &op_energy); 665*5754ecacSJeremy L Thompson CeedOperatorSetField(op_energy, "du", data[fine_level]->elem_restr_u, 666*5754ecacSJeremy L Thompson data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 667*5754ecacSJeremy L Thompson CeedOperatorSetField(op_energy, "geo_data", 668*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_i, 669*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, data[fine_level]->geo_data); 670*5754ecacSJeremy L Thompson CeedOperatorSetField(op_energy, "energy", data[fine_level]->elem_restr_energy, 671*5754ecacSJeremy L Thompson data[fine_level]->basis_energy, CEED_VECTOR_ACTIVE); 672*5754ecacSJeremy L Thompson // -- Save libCEED data 673*5754ecacSJeremy L Thompson data[fine_level]->qf_energy = qf_energy; 674*5754ecacSJeremy L Thompson data[fine_level]->op_energy = op_energy; 675*5754ecacSJeremy L Thompson 676*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 677*5754ecacSJeremy L Thompson // Diagnostic value computation 678*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 679*5754ecacSJeremy L Thompson // Create the QFunction and Operator that computes nodal diagnostic quantities 680*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 681*5754ecacSJeremy L Thompson // Geometric factors 682*5754ecacSJeremy L Thompson // -- Coordinate basis 683*5754ecacSJeremy L Thompson CeedBasis basis_x; 684*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS_LOBATTO, 685*5754ecacSJeremy L Thompson &basis_x); 686*5754ecacSJeremy L Thompson // -- QFunction 687*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.setup_geo, 688*5754ecacSJeremy L Thompson problem_data.setup_geo_loc, &qf_setup_geo); 689*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*dim, CEED_EVAL_GRAD); 690*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT); 691*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_setup_geo, "geo_data", geo_data_size, CEED_EVAL_NONE); 692*5754ecacSJeremy L Thompson // -- Operator 693*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, 694*5754ecacSJeremy L Thompson CEED_QFUNCTION_NONE, &op_setup_geo); 695*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "dx", data[fine_level]->elem_restr_x, 696*5754ecacSJeremy L Thompson basis_x, CEED_VECTOR_ACTIVE); 697*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, 698*5754ecacSJeremy L Thompson basis_x, CEED_VECTOR_NONE); 699*5754ecacSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "geo_data", 700*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_diagnostic_i, 701*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 702*5754ecacSJeremy L Thompson // -- Compute the quadrature data 703*5754ecacSJeremy L Thompson CeedOperatorApply(op_setup_geo, x_coord, data[fine_level]->geo_data_diagnostic, 704*5754ecacSJeremy L Thompson CEED_REQUEST_IMMEDIATE); 705*5754ecacSJeremy L Thompson // -- Cleanup 706*5754ecacSJeremy L Thompson CeedBasisDestroy(&basis_x); 707*5754ecacSJeremy L Thompson CeedQFunctionDestroy(&qf_setup_geo); 708*5754ecacSJeremy L Thompson CeedOperatorDestroy(&op_setup_geo); 709*5754ecacSJeremy L Thompson 710*5754ecacSJeremy L Thompson // Diagnostic quantities 711*5754ecacSJeremy L Thompson // -- QFunction 712*5754ecacSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem_data.diagnostic, 713*5754ecacSJeremy L Thompson problem_data.diagnostic_loc, &qf_diagnostic); 714*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_diagnostic, "u", num_comp_u, CEED_EVAL_INTERP); 715*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_diagnostic, "du", num_comp_u*dim, CEED_EVAL_GRAD); 716*5754ecacSJeremy L Thompson CeedQFunctionAddInput(qf_diagnostic, "geo_data", geo_data_size, CEED_EVAL_NONE); 717*5754ecacSJeremy L Thompson CeedQFunctionAddOutput(qf_diagnostic, "diagnostic", num_comp_u + num_comp_d, 718*5754ecacSJeremy L Thompson CEED_EVAL_NONE); 719*5754ecacSJeremy L Thompson CeedQFunctionSetContext(qf_diagnostic, phys_ctx); 720*5754ecacSJeremy L Thompson // -- Operator 721*5754ecacSJeremy L Thompson CeedOperatorCreate(ceed, qf_diagnostic, CEED_QFUNCTION_NONE, 722*5754ecacSJeremy L Thompson CEED_QFUNCTION_NONE, &op_diagnostic); 723*5754ecacSJeremy L Thompson CeedOperatorSetField(op_diagnostic, "u", data[fine_level]->elem_restr_u, 724*5754ecacSJeremy L Thompson data[fine_level]->basis_diagnostic, CEED_VECTOR_ACTIVE); 725*5754ecacSJeremy L Thompson CeedOperatorSetField(op_diagnostic, "du", data[fine_level]->elem_restr_u, 726*5754ecacSJeremy L Thompson data[fine_level]->basis_diagnostic, CEED_VECTOR_ACTIVE); 727*5754ecacSJeremy L Thompson CeedOperatorSetField(op_diagnostic, "geo_data", 728*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_geo_data_diagnostic_i, 729*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, data[fine_level]->geo_data_diagnostic); 730*5754ecacSJeremy L Thompson CeedOperatorSetField(op_diagnostic, "diagnostic", 731*5754ecacSJeremy L Thompson data[fine_level]->elem_restr_diagnostic, 732*5754ecacSJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 733*5754ecacSJeremy L Thompson // -- Save libCEED data 734*5754ecacSJeremy L Thompson data[fine_level]->qf_diagnostic = qf_diagnostic; 735*5754ecacSJeremy L Thompson data[fine_level]->op_diagnostic = op_diagnostic; 736*5754ecacSJeremy L Thompson 737*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 738*5754ecacSJeremy L Thompson // Cleanup 739*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 740*5754ecacSJeremy L Thompson CeedVectorDestroy(&x_coord); 741*5754ecacSJeremy L Thompson 742*5754ecacSJeremy L Thompson PetscFunctionReturn(0); 743*5754ecacSJeremy L Thompson }; 744*5754ecacSJeremy L Thompson 745*5754ecacSJeremy L Thompson // Set up libCEED multigrid level for a given degree 746*5754ecacSJeremy L Thompson // Prolongation and Restriction are between level and level+1 747*5754ecacSJeremy L Thompson PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx app_ctx, 748*5754ecacSJeremy L Thompson ProblemData problem_data, PetscInt level, 749*5754ecacSJeremy L Thompson PetscInt num_comp_u, PetscInt U_g_size, 750*5754ecacSJeremy L Thompson PetscInt U_loc_size, CeedVector fine_mult, 751*5754ecacSJeremy L Thompson CeedData *data) { 752*5754ecacSJeremy L Thompson PetscErrorCode ierr; 753*5754ecacSJeremy L Thompson CeedInt fine_level = app_ctx->num_levels - 1; 754*5754ecacSJeremy L Thompson CeedInt P = app_ctx->level_degrees[level] + 1; 755*5754ecacSJeremy L Thompson CeedInt Q = app_ctx->level_degrees[fine_level] + 1 + app_ctx->q_extra; 756*5754ecacSJeremy L Thompson CeedInt dim; 757*5754ecacSJeremy L Thompson CeedOperator op_jacobian, op_prolong, op_restrict; 758*5754ecacSJeremy L Thompson 759*5754ecacSJeremy L Thompson PetscFunctionBeginUser; 760*5754ecacSJeremy L Thompson 761*5754ecacSJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 762*5754ecacSJeremy L Thompson 763*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 764*5754ecacSJeremy L Thompson // libCEED restrictions 765*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 766*5754ecacSJeremy L Thompson // -- Solution restriction 767*5754ecacSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, P, 0, 0, 0, 768*5754ecacSJeremy L Thompson &data[level]->elem_restr_u); 769*5754ecacSJeremy L Thompson CHKERRQ(ierr); 770*5754ecacSJeremy L Thompson 771*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 772*5754ecacSJeremy L Thompson // libCEED bases 773*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 774*5754ecacSJeremy L Thompson // -- Solution basis 775*5754ecacSJeremy L Thompson CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, 776*5754ecacSJeremy L Thompson problem_data.quadrature_mode, 777*5754ecacSJeremy L Thompson &data[level]->basis_u); 778*5754ecacSJeremy L Thompson 779*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 780*5754ecacSJeremy L Thompson // Persistent libCEED vectors 781*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 782*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, U_loc_size, &data[level]->x_ceed); 783*5754ecacSJeremy L Thompson CeedVectorCreate(ceed, U_loc_size, &data[level]->y_ceed); 784*5754ecacSJeremy L Thompson 785*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 786*5754ecacSJeremy L Thompson // Coarse Grid, Prolongation, and Restriction Operators 787*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 788*5754ecacSJeremy L Thompson // Create the Operators that compute the prolongation and 789*5754ecacSJeremy L Thompson // restriction between the p-multigrid levels and the coarse grid eval. 790*5754ecacSJeremy L Thompson // --------------------------------------------------------------------------- 791*5754ecacSJeremy L Thompson CeedOperatorMultigridLevelCreate(data[level+1]->op_jacobian, fine_mult, 792*5754ecacSJeremy L Thompson data[level]->elem_restr_u, data[level]->basis_u, 793*5754ecacSJeremy L Thompson &op_jacobian, &op_prolong, &op_restrict); 794*5754ecacSJeremy L Thompson 795*5754ecacSJeremy L Thompson // -- Save libCEED data 796*5754ecacSJeremy L Thompson data[level]->op_jacobian = op_jacobian; 797*5754ecacSJeremy L Thompson data[level+1]->op_prolong = op_prolong; 798*5754ecacSJeremy L Thompson data[level+1]->op_restrict = op_restrict; 799*5754ecacSJeremy L Thompson 800*5754ecacSJeremy L Thompson PetscFunctionReturn(0); 801*5754ecacSJeremy L Thompson }; 802