13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 33d8e8822SJeremy L Thompson // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53d8e8822SJeremy L Thompson // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73d8e8822SJeremy L Thompson 85754ecacSJeremy L Thompson /// @file 95754ecacSJeremy L Thompson /// Command line option processing for solid mechanics example using PETSc 105754ecacSJeremy L Thompson 115754ecacSJeremy L Thompson #include "../include/cl-options.h" 125754ecacSJeremy L Thompson 135754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 145754ecacSJeremy L Thompson // Process command line options 155754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 165754ecacSJeremy L Thompson // Process general command line options 175754ecacSJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx) { 185754ecacSJeremy L Thompson PetscBool ceed_flag = PETSC_FALSE; 195754ecacSJeremy L Thompson 205754ecacSJeremy L Thompson PetscFunctionBeginUser; 215754ecacSJeremy L Thompson 222b730f8bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Elasticity / Hyperelasticity in PETSc with libCEED", NULL); 235754ecacSJeremy L Thompson 242b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource, 252b730f8bSJeremy L Thompson sizeof(app_ctx->ceed_resource), &ceed_flag)); 265754ecacSJeremy L Thompson 272b730f8bSJeremy L Thompson PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2)); // Default - current directory 282b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL)); 295754ecacSJeremy L Thompson 305754ecacSJeremy L Thompson app_ctx->degree = 3; 312b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, app_ctx->degree, &app_ctx->degree, NULL)); 325754ecacSJeremy L Thompson 335754ecacSJeremy L Thompson app_ctx->q_extra = 0; 342b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL)); 355754ecacSJeremy L Thompson 362b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, app_ctx->mesh_file, app_ctx->mesh_file, sizeof(app_ctx->mesh_file), NULL)); 375754ecacSJeremy L Thompson 385754ecacSJeremy L Thompson app_ctx->problem_choice = ELAS_LINEAR; // Default - Linear Elasticity 392b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-problem", "Solves Elasticity & Hyperelasticity Problems", NULL, problemTypes, (PetscEnum)app_ctx->problem_choice, 402b730f8bSJeremy L Thompson (PetscEnum *)&app_ctx->problem_choice, NULL)); 415754ecacSJeremy L Thompson app_ctx->name = problemTypes[app_ctx->problem_choice]; 425754ecacSJeremy L Thompson app_ctx->name_for_disp = problemTypesForDisp[app_ctx->problem_choice]; 435754ecacSJeremy L Thompson 445754ecacSJeremy L Thompson app_ctx->num_increments = app_ctx->problem_choice == ELAS_LINEAR ? 1 : 10; 452b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-num_steps", "Number of pseudo-time steps", NULL, app_ctx->num_increments, &app_ctx->num_increments, NULL)); 465754ecacSJeremy L Thompson 475754ecacSJeremy L Thompson app_ctx->forcing_choice = FORCE_NONE; // Default - no forcing term 482b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-forcing", "Set forcing function option", NULL, forcing_types, (PetscEnum)app_ctx->forcing_choice, 492b730f8bSJeremy L Thompson (PetscEnum *)&app_ctx->forcing_choice, NULL)); 505754ecacSJeremy L Thompson 515754ecacSJeremy L Thompson PetscInt max_n = 3; 525754ecacSJeremy L Thompson app_ctx->forcing_vector[0] = 0; 535754ecacSJeremy L Thompson app_ctx->forcing_vector[1] = -1; 545754ecacSJeremy L Thompson app_ctx->forcing_vector[2] = 0; 552b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalarArray("-forcing_vec", "Direction to apply constant force", NULL, app_ctx->forcing_vector, &max_n, NULL)); 565754ecacSJeremy L Thompson 572b730f8bSJeremy L Thompson if ((app_ctx->problem_choice == ELAS_FSInitial_NH1 || app_ctx->problem_choice == ELAS_FSInitial_NH2 || 582b730f8bSJeremy L Thompson app_ctx->problem_choice == ELAS_FSCurrent_NH1 || app_ctx->problem_choice == ELAS_FSCurrent_NH2 || 595754ecacSJeremy L Thompson app_ctx->problem_choice == ELAS_FSInitial_MR1) && 605754ecacSJeremy L Thompson app_ctx->forcing_choice == FORCE_CONST) 615754ecacSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, 625754ecacSJeremy L Thompson "Cannot use constant forcing and finite strain formulation. " 635754ecacSJeremy L Thompson "Constant forcing in reference frame currently unavaliable."); 645754ecacSJeremy L Thompson 655754ecacSJeremy L Thompson // Dirichlet boundary conditions 665754ecacSJeremy L Thompson app_ctx->bc_clamp_count = 16; 672b730f8bSJeremy L Thompson PetscCall( 682b730f8bSJeremy L Thompson PetscOptionsIntArray("-bc_clamp", "Face IDs to apply incremental Dirichlet BC", NULL, app_ctx->bc_clamp_faces, &app_ctx->bc_clamp_count, NULL)); 695754ecacSJeremy L Thompson // Set vector for each clamped BC 705754ecacSJeremy L Thompson for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) { 715754ecacSJeremy L Thompson // Translation vector 725754ecacSJeremy L Thompson char option_name[25]; 732b730f8bSJeremy L Thompson const size_t nclamp_params = sizeof(app_ctx->bc_clamp_max[0]) / sizeof(app_ctx->bc_clamp_max[0][0]); 742b730f8bSJeremy L Thompson for (PetscInt j = 0; j < nclamp_params; j++) app_ctx->bc_clamp_max[i][j] = 0.; 755754ecacSJeremy L Thompson 762b730f8bSJeremy L Thompson snprintf(option_name, sizeof option_name, "-bc_clamp_%" PetscInt_FMT "_translate", app_ctx->bc_clamp_faces[i]); 775754ecacSJeremy L Thompson max_n = 3; 782b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalarArray(option_name, "Vector to translate clamped end by", NULL, app_ctx->bc_clamp_max[i], &max_n, NULL)); 795754ecacSJeremy L Thompson 805754ecacSJeremy L Thompson // Rotation vector 815754ecacSJeremy L Thompson max_n = 5; 822b730f8bSJeremy L Thompson snprintf(option_name, sizeof option_name, "-bc_clamp_%" PetscInt_FMT "_rotate", app_ctx->bc_clamp_faces[i]); 832b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalarArray(option_name, "Vector with axis of rotation and rotation, in radians", NULL, &app_ctx->bc_clamp_max[i][3], 842b730f8bSJeremy L Thompson &max_n, NULL)); 855754ecacSJeremy L Thompson 865754ecacSJeremy L Thompson // Normalize 872b730f8bSJeremy L Thompson PetscScalar norm = sqrt(app_ctx->bc_clamp_max[i][3] * app_ctx->bc_clamp_max[i][3] + app_ctx->bc_clamp_max[i][4] * app_ctx->bc_clamp_max[i][4] + 882b730f8bSJeremy L Thompson app_ctx->bc_clamp_max[i][5] * app_ctx->bc_clamp_max[i][5]); 892b730f8bSJeremy L Thompson if (fabs(norm) < 1e-16) norm = 1; 902b730f8bSJeremy L Thompson for (PetscInt j = 0; j < 3; j++) app_ctx->bc_clamp_max[i][3 + j] /= norm; 915754ecacSJeremy L Thompson } 925754ecacSJeremy L Thompson 935754ecacSJeremy L Thompson // Neumann boundary conditions 945754ecacSJeremy L Thompson app_ctx->bc_traction_count = 16; 952b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_traction", "Face IDs to apply traction (Neumann) BC", NULL, app_ctx->bc_traction_faces, 962b730f8bSJeremy L Thompson &app_ctx->bc_traction_count, NULL)); 975754ecacSJeremy L Thompson // Set vector for each traction BC 985754ecacSJeremy L Thompson for (PetscInt i = 0; i < app_ctx->bc_traction_count; i++) { 995754ecacSJeremy L Thompson // Translation vector 1005754ecacSJeremy L Thompson char option_name[25]; 1012b730f8bSJeremy L Thompson for (PetscInt j = 0; j < 3; j++) app_ctx->bc_traction_vector[i][j] = 0.; 1025754ecacSJeremy L Thompson 1032b730f8bSJeremy L Thompson snprintf(option_name, sizeof option_name, "-bc_traction_%" PetscInt_FMT, app_ctx->bc_traction_faces[i]); 1045754ecacSJeremy L Thompson max_n = 3; 1055754ecacSJeremy L Thompson PetscBool set = false; 1062b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalarArray(option_name, "Traction vector for constrained face", NULL, app_ctx->bc_traction_vector[i], &max_n, &set)); 1075754ecacSJeremy L Thompson 1082b730f8bSJeremy L Thompson if (!set) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Traction vector must be set for all traction boundary conditions."); 1095754ecacSJeremy L Thompson } 1105754ecacSJeremy L Thompson 1115754ecacSJeremy L Thompson app_ctx->multigrid_choice = MULTIGRID_LOGARITHMIC; 1122b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-multigrid", "Set multigrid type option", NULL, multigrid_types, (PetscEnum)app_ctx->multigrid_choice, 1132b730f8bSJeremy L Thompson (PetscEnum *)&app_ctx->multigrid_choice, NULL)); 1145754ecacSJeremy L Thompson 1155754ecacSJeremy L Thompson app_ctx->test_mode = PETSC_FALSE; 1162b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, app_ctx->test_mode, &(app_ctx->test_mode), NULL)); 1175754ecacSJeremy L Thompson 1185754ecacSJeremy L Thompson app_ctx->expect_final_strain = -1.; 1192b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-expect_final_strain_energy", "Expect final strain energy close to this value.", NULL, app_ctx->expect_final_strain, 1202b730f8bSJeremy L Thompson &app_ctx->expect_final_strain, NULL)); 1215754ecacSJeremy L Thompson 1225754ecacSJeremy L Thompson app_ctx->test_tol = 1e-8; 1232b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-expect_final_state_rtol", "Relative tolerance for final strain energy test", NULL, app_ctx->test_tol, 1242b730f8bSJeremy L Thompson &app_ctx->test_tol, NULL)); 1255754ecacSJeremy L Thompson 1265754ecacSJeremy L Thompson app_ctx->view_soln = PETSC_FALSE; 1272b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-view_soln", "Write out solution vector for viewing", NULL, app_ctx->view_soln, &(app_ctx->view_soln), NULL)); 1285754ecacSJeremy L Thompson 1295754ecacSJeremy L Thompson app_ctx->view_final_soln = PETSC_FALSE; 1302b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-view_final_soln", "Write out final solution vector for viewing", NULL, app_ctx->view_final_soln, 1312b730f8bSJeremy L Thompson &(app_ctx->view_final_soln), NULL)); 1325754ecacSJeremy L Thompson 1335754ecacSJeremy L Thompson PetscBool set; 1345754ecacSJeremy L Thompson char energy_viewer_filename[PETSC_MAX_PATH_LEN] = ""; 1352b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-strain_energy_monitor", "Print out current strain energy at every load increment", NULL, energy_viewer_filename, 1362b730f8bSJeremy L Thompson energy_viewer_filename, sizeof(energy_viewer_filename), &set)); 1375754ecacSJeremy L Thompson if (set) { 1382b730f8bSJeremy L Thompson PetscCall(PetscViewerASCIIOpen(comm, energy_viewer_filename, &app_ctx->energy_viewer)); 1392b730f8bSJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "increment,energy\n")); 140*ea61e9acSJeremy L Thompson // Initial configuration is base energy state; this may not be true if we extend in the future to initially loaded configurations (because a truly 141*ea61e9acSJeremy L Thompson // at-rest initial state may not be realizable). 1422b730f8bSJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", 0., 0.)); 1435754ecacSJeremy L Thompson } 14467490bc6SJeremy L Thompson PetscOptionsEnd(); // End of setting AppCtx 1455754ecacSJeremy L Thompson 1465754ecacSJeremy L Thompson // Check for all required values set 1475754ecacSJeremy L Thompson if (app_ctx->test_mode) { 1482b730f8bSJeremy L Thompson if (app_ctx->forcing_choice == FORCE_NONE && !app_ctx->bc_clamp_count) app_ctx->forcing_choice = FORCE_MMS; 1495754ecacSJeremy L Thompson } 1505754ecacSJeremy L Thompson if (!app_ctx->bc_clamp_count && app_ctx->forcing_choice != FORCE_MMS) { 1515754ecacSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "-boundary options needed"); 1525754ecacSJeremy L Thompson } 1535754ecacSJeremy L Thompson 1545754ecacSJeremy L Thompson // Provide default ceed resource if not specified 1555754ecacSJeremy L Thompson if (!ceed_flag) { 1565754ecacSJeremy L Thompson const char *ceed_resource = "/cpu/self"; 1575754ecacSJeremy L Thompson strncpy(app_ctx->ceed_resource, ceed_resource, 10); 1585754ecacSJeremy L Thompson } 1595754ecacSJeremy L Thompson 1605754ecacSJeremy L Thompson // Determine number of levels 1615754ecacSJeremy L Thompson switch (app_ctx->multigrid_choice) { 1625754ecacSJeremy L Thompson case MULTIGRID_LOGARITHMIC: 1635754ecacSJeremy L Thompson app_ctx->num_levels = ceil(log(app_ctx->degree) / log(2)) + 1; 1645754ecacSJeremy L Thompson break; 1655754ecacSJeremy L Thompson case MULTIGRID_UNIFORM: 1665754ecacSJeremy L Thompson app_ctx->num_levels = app_ctx->degree; 1675754ecacSJeremy L Thompson break; 1685754ecacSJeremy L Thompson case MULTIGRID_NONE: 1695754ecacSJeremy L Thompson app_ctx->num_levels = 1; 1705754ecacSJeremy L Thompson break; 1715754ecacSJeremy L Thompson } 1725754ecacSJeremy L Thompson 1735754ecacSJeremy L Thompson // Populate array of degrees for each level for multigrid 1742b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(app_ctx->num_levels, &(app_ctx->level_degrees))); 1755754ecacSJeremy L Thompson 1765754ecacSJeremy L Thompson switch (app_ctx->multigrid_choice) { 1775754ecacSJeremy L Thompson case MULTIGRID_LOGARITHMIC: 1782b730f8bSJeremy L Thompson for (int i = 0; i < app_ctx->num_levels - 1; i++) app_ctx->level_degrees[i] = pow(2, i); 1795754ecacSJeremy L Thompson app_ctx->level_degrees[app_ctx->num_levels - 1] = app_ctx->degree; 1805754ecacSJeremy L Thompson break; 1815754ecacSJeremy L Thompson case MULTIGRID_UNIFORM: 1822b730f8bSJeremy L Thompson for (int i = 0; i < app_ctx->num_levels; i++) app_ctx->level_degrees[i] = i + 1; 1835754ecacSJeremy L Thompson break; 1845754ecacSJeremy L Thompson case MULTIGRID_NONE: 1855754ecacSJeremy L Thompson app_ctx->level_degrees[0] = app_ctx->degree; 1865754ecacSJeremy L Thompson break; 1875754ecacSJeremy L Thompson } 1885754ecacSJeremy L Thompson 1895754ecacSJeremy L Thompson PetscFunctionReturn(0); 1905754ecacSJeremy L Thompson }; 191