xref: /libCEED/examples/solids/src/matops.c (revision cd11335e55ed5306fc93608d50d5700c3a391b92)
1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4ccaff030SJeremy L Thompson //
5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson /// @file
18ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc
19ccaff030SJeremy L Thompson 
20ccaff030SJeremy L Thompson #include "../elasticity.h"
21ccaff030SJeremy L Thompson 
22ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
23ccaff030SJeremy L Thompson // libCEED Operators for MatShell
24ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
25ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
26ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) {
27ccaff030SJeremy L Thompson   PetscErrorCode ierr;
28ccaff030SJeremy L Thompson   PetscScalar *x, *y;
29ccaff030SJeremy L Thompson 
30ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
31ccaff030SJeremy L Thompson 
32ccaff030SJeremy L Thompson   // Global-to-local
33ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
34ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Yloc); CHKERRQ(ierr);
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson   // Setup CEED vectors
3762e9c006SJeremy L Thompson   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
3862e9c006SJeremy L Thompson   CHKERRQ(ierr);
3962e9c006SJeremy L Thompson   ierr = user->VecGetArray(user->Yloc, &y); CHKERRQ(ierr);
4062e9c006SJeremy L Thompson   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
4162e9c006SJeremy L Thompson   CeedVectorSetArray(user->Yceed, user->memType, CEED_USE_POINTER, y);
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson   // Apply CEED operator
44ccaff030SJeremy L Thompson   CeedOperatorApply(user->op, user->Xceed, user->Yceed, CEED_REQUEST_IMMEDIATE);
4562e9c006SJeremy L Thompson   CeedVectorSyncArray(user->Yceed, user->memType);
46ccaff030SJeremy L Thompson 
47ccaff030SJeremy L Thompson   // Restore PETSc vectors
4862e9c006SJeremy L Thompson   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
49ccaff030SJeremy L Thompson   CHKERRQ(ierr);
5062e9c006SJeremy L Thompson   ierr = user->VecRestoreArray(user->Yloc, &y); CHKERRQ(ierr);
51ccaff030SJeremy L Thompson 
52ccaff030SJeremy L Thompson   // Local-to-global
53ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
54ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Yloc, ADD_VALUES, Y); CHKERRQ(ierr);
55ccaff030SJeremy L Thompson 
56ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
57ccaff030SJeremy L Thompson };
58ccaff030SJeremy L Thompson 
59ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
60ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
61ccaff030SJeremy L Thompson   PetscErrorCode ierr;
62ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
63ccaff030SJeremy L Thompson 
64ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
65ccaff030SJeremy L Thompson 
66ccaff030SJeremy L Thompson   // Use computed BCs
67ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
68ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
69ccaff030SJeremy L Thompson                                     user->loadIncrement, NULL, NULL, NULL);
70ccaff030SJeremy L Thompson   CHKERRQ(ierr);
71ccaff030SJeremy L Thompson 
72ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
73ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
74ccaff030SJeremy L Thompson 
75ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
76ccaff030SJeremy L Thompson };
77ccaff030SJeremy L Thompson 
78ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
79ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) {
80ccaff030SJeremy L Thompson   PetscErrorCode ierr;
81ccaff030SJeremy L Thompson   UserMult user = (UserMult)ctx;
82ccaff030SJeremy L Thompson 
83ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
84ccaff030SJeremy L Thompson 
8562e9c006SJeremy L Thompson   // Zero boundary values
86ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
87ccaff030SJeremy L Thompson 
88ccaff030SJeremy L Thompson   // libCEED for local action of residual evaluator
89ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
90ccaff030SJeremy L Thompson 
91ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
92ccaff030SJeremy L Thompson };
93ccaff030SJeremy L Thompson 
94ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
95ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) {
96ccaff030SJeremy L Thompson   PetscErrorCode ierr;
97ccaff030SJeremy L Thompson   UserMult user;
98ccaff030SJeremy L Thompson 
99ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
100ccaff030SJeremy L Thompson 
101ccaff030SJeremy L Thompson   // Zero boundary values
102ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
103ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
104ccaff030SJeremy L Thompson 
105ccaff030SJeremy L Thompson   // libCEED for local action of Jacobian
106ccaff030SJeremy L Thompson   ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr);
107ccaff030SJeremy L Thompson 
108ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
109ccaff030SJeremy L Thompson };
110ccaff030SJeremy L Thompson 
111ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
112ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) {
113ccaff030SJeremy L Thompson   PetscErrorCode ierr;
114ccaff030SJeremy L Thompson   UserMultProlongRestr user;
115ccaff030SJeremy L Thompson   PetscScalar *c, *f;
116ccaff030SJeremy L Thompson 
117ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
118ccaff030SJeremy L Thompson 
119ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
120ccaff030SJeremy L Thompson 
121ccaff030SJeremy L Thompson   // Global-to-local
122ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
123ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmC, X, INSERT_VALUES, user->locVecC);
124ccaff030SJeremy L Thompson   CHKERRQ(ierr);
125ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
126ccaff030SJeremy L Thompson 
127ccaff030SJeremy L Thompson   // Setup CEED vectors
12862e9c006SJeremy L Thompson   ierr = user->VecGetArrayRead(user->locVecC, (const PetscScalar **)&c);
129ccaff030SJeremy L Thompson   CHKERRQ(ierr);
13062e9c006SJeremy L Thompson   ierr = user->VecGetArray(user->locVecF, &f); CHKERRQ(ierr);
13162e9c006SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, user->memType, CEED_USE_POINTER, c);
13262e9c006SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, user->memType, CEED_USE_POINTER, f);
133ccaff030SJeremy L Thompson 
134ccaff030SJeremy L Thompson   // Apply CEED operator
135ccaff030SJeremy L Thompson   CeedOperatorApply(user->opProlong, user->ceedVecC, user->ceedVecF,
136ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
13762e9c006SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecF, user->memType);
138ccaff030SJeremy L Thompson 
139ccaff030SJeremy L Thompson   // Restore PETSc vectors
140*cd11335eSJeremy L Thompson   ierr = user->VecRestoreArrayRead(user->locVecC, (const PetscScalar **)&c);
141ccaff030SJeremy L Thompson   CHKERRQ(ierr);
14262e9c006SJeremy L Thompson   ierr = user->VecRestoreArray(user->locVecF, &f); CHKERRQ(ierr);
143ccaff030SJeremy L Thompson 
144ccaff030SJeremy L Thompson   // Multiplicity
145ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
146ccaff030SJeremy L Thompson 
147ccaff030SJeremy L Thompson   // Local-to-global
148ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
149ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmF, user->locVecF, ADD_VALUES, Y);
150ccaff030SJeremy L Thompson   CHKERRQ(ierr);
151ccaff030SJeremy L Thompson 
152ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
153ccaff030SJeremy L Thompson }
154ccaff030SJeremy L Thompson 
155ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
156ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) {
157ccaff030SJeremy L Thompson   PetscErrorCode ierr;
158ccaff030SJeremy L Thompson   UserMultProlongRestr user;
159ccaff030SJeremy L Thompson   PetscScalar *c, *f;
160ccaff030SJeremy L Thompson 
161ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
162ccaff030SJeremy L Thompson 
163ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
164ccaff030SJeremy L Thompson 
165ccaff030SJeremy L Thompson   // Global-to-local
166ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecF); CHKERRQ(ierr);
167ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dmF, X, INSERT_VALUES, user->locVecF);
168ccaff030SJeremy L Thompson   CHKERRQ(ierr);
169ccaff030SJeremy L Thompson   ierr = VecZeroEntries(user->locVecC); CHKERRQ(ierr);
170ccaff030SJeremy L Thompson 
171ccaff030SJeremy L Thompson   // Multiplicity
172ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(user->locVecF, user->locVecF, user->multVec);
173ccaff030SJeremy L Thompson   CHKERRQ(ierr);
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson   // Setup CEED vectors
17662e9c006SJeremy L Thompson   ierr = user->VecGetArrayRead(user->locVecF, (const PetscScalar **)&f);
17762e9c006SJeremy L Thompson   CHKERRQ(ierr);
17862e9c006SJeremy L Thompson   ierr = user->VecGetArray(user->locVecC, &c); CHKERRQ(ierr);
17962e9c006SJeremy L Thompson   CeedVectorSetArray(user->ceedVecF, user->memType, CEED_USE_POINTER, f);
18062e9c006SJeremy L Thompson   CeedVectorSetArray(user->ceedVecC, user->memType, CEED_USE_POINTER, c);
181ccaff030SJeremy L Thompson 
182ccaff030SJeremy L Thompson   // Apply CEED operator
183ccaff030SJeremy L Thompson   CeedOperatorApply(user->opRestrict, user->ceedVecF, user->ceedVecC,
184ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
185ccaff030SJeremy L Thompson   CeedVectorSyncArray(user->ceedVecC, CEED_MEM_HOST);
186ccaff030SJeremy L Thompson 
187ccaff030SJeremy L Thompson   // Restore PETSc vectors
18862e9c006SJeremy L Thompson   ierr = user->VecRestoreArrayRead(user->locVecF, (const PetscScalar **)&f);
189ccaff030SJeremy L Thompson   CHKERRQ(ierr);
19062e9c006SJeremy L Thompson   ierr = user->VecRestoreArray(user->locVecC, &c); CHKERRQ(ierr);
191ccaff030SJeremy L Thompson 
192ccaff030SJeremy L Thompson   // Local-to-global
193ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
194ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dmC, user->locVecC, ADD_VALUES, Y);
195ccaff030SJeremy L Thompson   CHKERRQ(ierr);
196ccaff030SJeremy L Thompson 
197ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
198ccaff030SJeremy L Thompson };
1992d93065eSjeremylt 
200ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
201ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) {
202ccaff030SJeremy L Thompson   PetscErrorCode ierr;
203ccaff030SJeremy L Thompson   UserMult user;
204ccaff030SJeremy L Thompson 
205ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
206ccaff030SJeremy L Thompson 
207ccaff030SJeremy L Thompson   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
208ccaff030SJeremy L Thompson 
209f7b4142eSJeremy L Thompson   // -- Set physics context
210f7b4142eSJeremy L Thompson   if (user->physSmoother)
211f7b4142eSJeremy L Thompson     CeedQFunctionSetContext(user->qf, user->physSmoother,
212f7b4142eSJeremy L Thompson                             sizeof(*user->physSmoother));
213f7b4142eSJeremy L Thompson 
2142bba3ffaSJeremy L Thompson   // Compute Diagonal via libCEED
2152bba3ffaSJeremy L Thompson   PetscScalar *x;
2162bba3ffaSJeremy L Thompson 
2172bba3ffaSJeremy L Thompson   // -- Place PETSc vector in libCEED vector
21862e9c006SJeremy L Thompson   ierr = user->VecGetArray(user->Xloc, &x); CHKERRQ(ierr);
21962e9c006SJeremy L Thompson   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
2202bba3ffaSJeremy L Thompson 
2212bba3ffaSJeremy L Thompson   // -- Compute Diagonal
2222bba3ffaSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op, user->Xceed,
223ccaff030SJeremy L Thompson                                      CEED_REQUEST_IMMEDIATE);
22462e9c006SJeremy L Thompson   CeedVectorSyncArray(user->Xceed, user->memType);
225ccaff030SJeremy L Thompson 
226f7b4142eSJeremy L Thompson   // -- Reset physics context
227f7b4142eSJeremy L Thompson   if (user->physSmoother)
228f7b4142eSJeremy L Thompson     CeedQFunctionSetContext(user->qf, user->phys, sizeof(*user->phys));
229f7b4142eSJeremy L Thompson 
230ccaff030SJeremy L Thompson   // -- Local-to-Global
23162e9c006SJeremy L Thompson   ierr = user->VecRestoreArray(user->Xloc, &x); CHKERRQ(ierr);
232ccaff030SJeremy L Thompson   ierr = VecZeroEntries(D); CHKERRQ(ierr);
233ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, user->Xloc, ADD_VALUES, D); CHKERRQ(ierr);
234ccaff030SJeremy L Thompson 
2352bba3ffaSJeremy L Thompson   // Cleanup
2362bba3ffaSJeremy L Thompson   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
237ccaff030SJeremy L Thompson 
238ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
239ccaff030SJeremy L Thompson };
2402d93065eSjeremylt 
2412d93065eSjeremylt // This function calculates the strain energy in the final solution
242a3c02c40SJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
243a3c02c40SJeremy L Thompson                                    CeedOperator opEnergy, Vec X,
244a3c02c40SJeremy L Thompson                                    PetscReal *energy) {
2452d93065eSjeremylt   PetscErrorCode ierr;
2462d93065eSjeremylt   PetscScalar *x;
2472d93065eSjeremylt   CeedInt length;
2482d93065eSjeremylt 
2492d93065eSjeremylt   PetscFunctionBeginUser;
2502d93065eSjeremylt 
2512d93065eSjeremylt   // Global-to-local
2522d93065eSjeremylt   ierr = VecZeroEntries(user->Xloc); CHKERRQ(ierr);
2532d93065eSjeremylt   ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->Xloc); CHKERRQ(ierr);
2542d93065eSjeremylt   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->Xloc,
2552d93065eSjeremylt                                     user->loadIncrement, NULL, NULL, NULL);
2565c25879aSJeremy L Thompson   CHKERRQ(ierr);
2572d93065eSjeremylt 
258a3c02c40SJeremy L Thompson   // Setup libCEED input vector
25962e9c006SJeremy L Thompson   ierr = user->VecGetArrayRead(user->Xloc, (const PetscScalar **)&x);
2602d93065eSjeremylt   CHKERRQ(ierr);
26162e9c006SJeremy L Thompson   CeedVectorSetArray(user->Xceed, user->memType, CEED_USE_POINTER, x);
2622d93065eSjeremylt 
263a3c02c40SJeremy L Thompson   // Setup libCEED output vector
264a3c02c40SJeremy L Thompson   Vec Eloc;
265a3c02c40SJeremy L Thompson   CeedVector eloc;
266a3c02c40SJeremy L Thompson   ierr = DMCreateLocalVector(dmEnergy, &Eloc); CHKERRQ(ierr);
267a3c02c40SJeremy L Thompson   ierr = VecGetSize(Eloc, &length); CHKERRQ(ierr);
268a3c02c40SJeremy L Thompson   ierr = VecDestroy(&Eloc); CHKERRQ(ierr);
269a3c02c40SJeremy L Thompson   CeedVectorCreate(user->ceed, length, &eloc);
270a3c02c40SJeremy L Thompson 
2712d93065eSjeremylt   // Apply libCEED operator
272a3c02c40SJeremy L Thompson   CeedOperatorApply(opEnergy, user->Xceed, eloc, CEED_REQUEST_IMMEDIATE);
2732d93065eSjeremylt 
2742d93065eSjeremylt   // Restore PETSc vector
27562e9c006SJeremy L Thompson   ierr = user->VecRestoreArrayRead(user->Xloc, (const PetscScalar **)&x);
2762d93065eSjeremylt   CHKERRQ(ierr);
2772d93065eSjeremylt 
2782d93065eSjeremylt   // Reduce max error
2792d93065eSjeremylt   const CeedScalar *e;
280a3c02c40SJeremy L Thompson   CeedVectorGetArrayRead(eloc, CEED_MEM_HOST, &e);
2812d93065eSjeremylt   (*energy) = 0;
2822d93065eSjeremylt   for (CeedInt i=0; i<length; i++)
2832d93065eSjeremylt     (*energy) += e[i];
284a3c02c40SJeremy L Thompson   CeedVectorRestoreArrayRead(eloc, &e);
285a3c02c40SJeremy L Thompson   CeedVectorDestroy(&eloc);
2862d93065eSjeremylt 
2872d93065eSjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM,
2882d93065eSjeremylt                        user->comm); CHKERRQ(ierr);
2892d93065eSjeremylt 
2902d93065eSjeremylt   PetscFunctionReturn(0);
2912d93065eSjeremylt };
292