xref: /libCEED/examples/solids/src/misc.c (revision ccaff0309dc399f656ea11018b919b30feb8b0fa)
1*ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2*ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3*ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4*ccaff030SJeremy L Thompson //
5*ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6*ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7*ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8*ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9*ccaff030SJeremy L Thompson //
10*ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11*ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12*ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13*ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14*ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15*ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16*ccaff030SJeremy L Thompson 
17*ccaff030SJeremy L Thompson /// @file
18*ccaff030SJeremy L Thompson /// Helper functions for solid mechanics example using PETSc
19*ccaff030SJeremy L Thompson 
20*ccaff030SJeremy L Thompson #include "../elasticity.h"
21*ccaff030SJeremy L Thompson 
22*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
23*ccaff030SJeremy L Thompson // Create libCEED operator context
24*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
25*ccaff030SJeremy L Thompson // Setup context data for Jacobian evaluation
26*ccaff030SJeremy L Thompson PetscErrorCode SetupJacobianCtx(MPI_Comm comm, AppCtx appCtx, DM dm, Vec V,
27*ccaff030SJeremy L Thompson                                 Vec Vloc, CeedData ceedData, Ceed ceed,
28*ccaff030SJeremy L Thompson                                 UserMult jacobianCtx) {
29*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
30*ccaff030SJeremy L Thompson 
31*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
32*ccaff030SJeremy L Thompson 
33*ccaff030SJeremy L Thompson   // PETSc objects
34*ccaff030SJeremy L Thompson   jacobianCtx->comm = comm;
35*ccaff030SJeremy L Thompson   jacobianCtx->dm = dm;
36*ccaff030SJeremy L Thompson 
37*ccaff030SJeremy L Thompson   // Work vectors
38*ccaff030SJeremy L Thompson   jacobianCtx->Xloc = Vloc;
39*ccaff030SJeremy L Thompson   ierr = VecDuplicate(Vloc, &jacobianCtx->Yloc); CHKERRQ(ierr);
40*ccaff030SJeremy L Thompson   jacobianCtx->Xceed = ceedData->xceed;
41*ccaff030SJeremy L Thompson   jacobianCtx->Yceed = ceedData->yceed;
42*ccaff030SJeremy L Thompson 
43*ccaff030SJeremy L Thompson   // libCEED operator
44*ccaff030SJeremy L Thompson   jacobianCtx->op = ceedData->opJacob;
45*ccaff030SJeremy L Thompson 
46*ccaff030SJeremy L Thompson   // Ceed
47*ccaff030SJeremy L Thompson   jacobianCtx->ceed = ceed;
48*ccaff030SJeremy L Thompson 
49*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
50*ccaff030SJeremy L Thompson };
51*ccaff030SJeremy L Thompson 
52*ccaff030SJeremy L Thompson // Setup context data for prolongation and restriction operators
53*ccaff030SJeremy L Thompson PetscErrorCode SetupProlongRestrictCtx(MPI_Comm comm, DM dmC, DM dmF, Vec VF,
54*ccaff030SJeremy L Thompson                                        Vec VlocC, Vec VlocF, CeedData ceedDataC,
55*ccaff030SJeremy L Thompson                                        CeedData ceedDataF, Ceed ceed,
56*ccaff030SJeremy L Thompson                                        UserMultProlongRestr prolongRestrCtx) {
57*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
58*ccaff030SJeremy L Thompson   PetscScalar *multArray;
59*ccaff030SJeremy L Thompson 
60*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
61*ccaff030SJeremy L Thompson 
62*ccaff030SJeremy L Thompson   // PETSc objects
63*ccaff030SJeremy L Thompson   prolongRestrCtx->comm = comm;
64*ccaff030SJeremy L Thompson   prolongRestrCtx->dmC = dmC;
65*ccaff030SJeremy L Thompson   prolongRestrCtx->dmF = dmF;
66*ccaff030SJeremy L Thompson 
67*ccaff030SJeremy L Thompson   // Work vectors
68*ccaff030SJeremy L Thompson   prolongRestrCtx->locVecC = VlocC;
69*ccaff030SJeremy L Thompson   prolongRestrCtx->locVecF = VlocF;
70*ccaff030SJeremy L Thompson   prolongRestrCtx->ceedVecC = ceedDataC->xceed;
71*ccaff030SJeremy L Thompson   prolongRestrCtx->ceedVecF = ceedDataF->xceed;
72*ccaff030SJeremy L Thompson 
73*ccaff030SJeremy L Thompson   // libCEED operators
74*ccaff030SJeremy L Thompson   prolongRestrCtx->opProlong = ceedDataF->opProlong;
75*ccaff030SJeremy L Thompson   prolongRestrCtx->opRestrict = ceedDataF->opRestrict;
76*ccaff030SJeremy L Thompson 
77*ccaff030SJeremy L Thompson   // Ceed
78*ccaff030SJeremy L Thompson   prolongRestrCtx->ceed = ceed;
79*ccaff030SJeremy L Thompson 
80*ccaff030SJeremy L Thompson   // Multiplicity vector
81*ccaff030SJeremy L Thompson   // -- Set libCEED vector
82*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(VlocF);
83*ccaff030SJeremy L Thompson   ierr = VecGetArray(VlocF, &multArray); CHKERRQ(ierr);
84*ccaff030SJeremy L Thompson   CeedVectorSetArray(ceedDataF->xceed, CEED_MEM_HOST, CEED_USE_POINTER,
85*ccaff030SJeremy L Thompson                      multArray);
86*ccaff030SJeremy L Thompson 
87*ccaff030SJeremy L Thompson   // -- Compute multiplicity
88*ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(ceedDataF->Erestrictu, ceedDataF->xceed);
89*ccaff030SJeremy L Thompson 
90*ccaff030SJeremy L Thompson   // -- Restore PETSc vector
91*ccaff030SJeremy L Thompson   ierr = VecRestoreArray(VlocF, &multArray); CHKERRQ(ierr);
92*ccaff030SJeremy L Thompson 
93*ccaff030SJeremy L Thompson   // -- Local-to-global
94*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(VF); CHKERRQ(ierr);
95*ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dmF, VlocF, ADD_VALUES, VF); CHKERRQ(ierr);
96*ccaff030SJeremy L Thompson 
97*ccaff030SJeremy L Thompson   // -- Global-to-local
98*ccaff030SJeremy L Thompson   ierr = VecDuplicate(VlocF, &prolongRestrCtx->multVec); CHKERRQ(ierr);
99*ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(dmF, VF, INSERT_VALUES, prolongRestrCtx->multVec);
100*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
101*ccaff030SJeremy L Thompson 
102*ccaff030SJeremy L Thompson   // -- Reciprocal
103*ccaff030SJeremy L Thompson   ierr = VecReciprocal(prolongRestrCtx->multVec); CHKERRQ(ierr);
104*ccaff030SJeremy L Thompson 
105*ccaff030SJeremy L Thompson   // -- Reset work arrays
106*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(VF); CHKERRQ(ierr);
107*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(VlocF); CHKERRQ(ierr);
108*ccaff030SJeremy L Thompson 
109*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
110*ccaff030SJeremy L Thompson };
111*ccaff030SJeremy L Thompson 
112*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
113*ccaff030SJeremy L Thompson // Jacobian setup
114*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
115*ccaff030SJeremy L Thompson PetscErrorCode FormJacobian(SNES snes, Vec U, Mat J, Mat Jpre, void *ctx) {
116*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
117*ccaff030SJeremy L Thompson 
118*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
119*ccaff030SJeremy L Thompson 
120*ccaff030SJeremy L Thompson   // Context data
121*ccaff030SJeremy L Thompson   FormJacobCtx  formJacobCtx = (FormJacobCtx)ctx;
122*ccaff030SJeremy L Thompson   PetscInt      numLevels = formJacobCtx->numLevels;
123*ccaff030SJeremy L Thompson   Mat           *jacobMat = formJacobCtx->jacobMat;
124*ccaff030SJeremy L Thompson 
125*ccaff030SJeremy L Thompson   // Update Jacobian PC smoother on each level
126*ccaff030SJeremy L Thompson   for (int level = 0; level < numLevels; level++) {
127*ccaff030SJeremy L Thompson     // -- Update diagonal state counter
128*ccaff030SJeremy L Thompson     ierr = MatAssemblyBegin(jacobMat[level], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
129*ccaff030SJeremy L Thompson     ierr = MatAssemblyEnd(jacobMat[level], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
130*ccaff030SJeremy L Thompson   }
131*ccaff030SJeremy L Thompson 
132*ccaff030SJeremy L Thompson   // Form coarse assembled matrix
133*ccaff030SJeremy L Thompson   ierr = VecZeroEntries(formJacobCtx->Ucoarse); CHKERRQ(ierr);
134*ccaff030SJeremy L Thompson   ierr = SNESComputeJacobianDefaultColor(formJacobCtx->snesCoarse,
135*ccaff030SJeremy L Thompson                                          formJacobCtx->Ucoarse,
136*ccaff030SJeremy L Thompson                                          formJacobCtx->jacobMat[0],
137*ccaff030SJeremy L Thompson                                          formJacobCtx->jacobMatCoarse, NULL);
138*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
139*ccaff030SJeremy L Thompson 
140*ccaff030SJeremy L Thompson   // Jpre might be AIJ (e.g., when using coloring), so we need to assemble it
141*ccaff030SJeremy L Thompson   ierr = MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
142*ccaff030SJeremy L Thompson   ierr = MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
143*ccaff030SJeremy L Thompson 
144*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
145*ccaff030SJeremy L Thompson };
146*ccaff030SJeremy L Thompson 
147*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
148*ccaff030SJeremy L Thompson // Output solution for visualization
149*ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
150*ccaff030SJeremy L Thompson PetscErrorCode ViewSolution(MPI_Comm comm, Vec U, PetscInt increment,
151*ccaff030SJeremy L Thompson                             PetscScalar loadIncrement) {
152*ccaff030SJeremy L Thompson   PetscErrorCode ierr;
153*ccaff030SJeremy L Thompson   DM dm;
154*ccaff030SJeremy L Thompson   PetscViewer viewer;
155*ccaff030SJeremy L Thompson   char outputFilename[PETSC_MAX_PATH_LEN];
156*ccaff030SJeremy L Thompson 
157*ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
158*ccaff030SJeremy L Thompson 
159*ccaff030SJeremy L Thompson   // Build file name
160*ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(outputFilename, sizeof outputFilename,
161*ccaff030SJeremy L Thompson                        "solution-%03D.vtu", increment); CHKERRQ(ierr);
162*ccaff030SJeremy L Thompson 
163*ccaff030SJeremy L Thompson   // Increment senquence
164*ccaff030SJeremy L Thompson   ierr = VecGetDM(U, &dm); CHKERRQ(ierr);
165*ccaff030SJeremy L Thompson   ierr = DMSetOutputSequenceNumber(dm, increment, loadIncrement); CHKERRQ(ierr);
166*ccaff030SJeremy L Thompson 
167*ccaff030SJeremy L Thompson   // Output solution vector
168*ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(comm, outputFilename, FILE_MODE_WRITE, &viewer);
169*ccaff030SJeremy L Thompson   CHKERRQ(ierr);
170*ccaff030SJeremy L Thompson   ierr = VecView(U, viewer); CHKERRQ(ierr);
171*ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
172*ccaff030SJeremy L Thompson 
173*ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
174*ccaff030SJeremy L Thompson };
175