xref: /libCEED/examples/solids/problems/neo-hookean.c (revision b2e3f8ecbfa285d0f4ffde9b24c57cc13f0319fb)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include "../problems/neo-hookean.h"
9 
10 #include <ceed.h>
11 #include <petscsys.h>
12 
13 // Build libCEED context object
14 PetscErrorCode PhysicsContext_NH(MPI_Comm comm, Ceed ceed, Units *units, CeedQFunctionContext *ctx) {
15   Physics_NH phys;
16 
17   PetscFunctionBegin;
18 
19   PetscCall(PetscMalloc1(1, units));
20   PetscCall(PetscMalloc1(1, &phys));
21   PetscCall(ProcessPhysics_NH(comm, phys, *units));
22   CeedQFunctionContextCreate(ceed, ctx);
23   CeedQFunctionContextSetData(*ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(*phys), phys);
24   PetscCall(PetscFree(phys));
25 
26   PetscFunctionReturn(0);
27 }
28 
29 // Build libCEED smoother context object
30 PetscErrorCode PhysicsSmootherContext_NH(MPI_Comm comm, Ceed ceed, CeedQFunctionContext ctx, CeedQFunctionContext *ctx_smoother) {
31   PetscScalar nu_smoother = 0;
32   PetscBool   nu_flag     = PETSC_FALSE;
33   Physics_NH  phys, phys_smoother;
34 
35   PetscFunctionBegin;
36 
37   PetscOptionsBegin(comm, NULL, "Neo-Hookean physical parameters for smoother", NULL);
38 
39   PetscCall(PetscOptionsScalar("-nu_smoother", "Poisson's ratio for smoother", NULL, nu_smoother, &nu_smoother, &nu_flag));
40 
41   PetscOptionsEnd();  // End of setting Physics
42 
43   if (nu_flag) {
44     // Copy context
45     CeedQFunctionContextGetData(ctx, CEED_MEM_HOST, &phys);
46     PetscCall(PetscMalloc1(1, &phys_smoother));
47     PetscCall(PetscMemcpy(phys_smoother, phys, sizeof(*phys)));
48     CeedQFunctionContextRestoreData(ctx, &phys);
49     // Create smoother context
50     CeedQFunctionContextCreate(ceed, ctx_smoother);
51     phys_smoother->nu = nu_smoother;
52     CeedQFunctionContextSetData(*ctx_smoother, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(*phys_smoother), phys_smoother);
53     PetscCall(PetscFree(phys_smoother));
54   } else {
55     *ctx_smoother = NULL;
56   }
57 
58   PetscFunctionReturn(0);
59 }
60 
61 // Process physics options - Neo-Hookean
62 PetscErrorCode ProcessPhysics_NH(MPI_Comm comm, Physics_NH phys, Units units) {
63   PetscBool nu_flag    = PETSC_FALSE;
64   PetscBool Young_flag = PETSC_FALSE;
65   phys->nu             = 0;
66   phys->E              = 0;
67   units->meter         = 1;  // 1 meter in scaled length units
68   units->second        = 1;  // 1 second in scaled time units
69   units->kilogram      = 1;  // 1 kilogram in scaled mass units
70 
71   PetscFunctionBeginUser;
72 
73   PetscOptionsBegin(comm, NULL, "Neo-Hookean physical parameters", NULL);
74 
75   PetscCall(PetscOptionsScalar("-nu", "Poisson's ratio", NULL, phys->nu, &phys->nu, &nu_flag));
76 
77   PetscCall(PetscOptionsScalar("-E", "Young's Modulus", NULL, phys->E, &phys->E, &Young_flag));
78 
79   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, units->meter, &units->meter, NULL));
80   units->meter = fabs(units->meter);
81 
82   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, units->second, &units->second, NULL));
83   units->second = fabs(units->second);
84 
85   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, units->kilogram, &units->kilogram, NULL));
86   units->kilogram = fabs(units->kilogram);
87 
88   PetscOptionsEnd();  // End of setting Physics
89 
90   // Check for all required options to be set
91   if (!nu_flag) {
92     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "-nu option needed");
93   }
94   if (!Young_flag) {
95     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "-E option needed");
96   }
97 
98   // Define derived units
99   units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
100 
101   // Scale E to Pa
102   phys->E *= units->Pascal;
103 
104   PetscFunctionReturn(0);
105 };