Actual source code: tsfd.c
2: #include <private/tsimpl.h> /*I "petscts.h" I*/
6: /*@C
7: TSDefaultComputeJacobianColor - Computes the Jacobian using
8: finite differences and coloring to exploit matrix sparsity.
9:
10: Collective on TS, Vec and Mat
12: Input Parameters:
13: + ts - nonlinear solver object
14: . t - current time
15: . x1 - location at which to evaluate Jacobian
16: - ctx - coloring context, where ctx must have type MatFDColoring,
17: as created via MatFDColoringCreate()
19: Output Parameters:
20: + J - Jacobian matrix (not altered in this routine)
21: . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
22: - flag - flag indicating whether the matrix sparsity structure has changed
24: Level: intermediate
26: .keywords: TS, finite differences, Jacobian, coloring, sparse
28: .seealso: TSSetJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
29: @*/
30: PetscErrorCode TSDefaultComputeJacobianColor(TS ts,PetscReal t,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
31: {
32: MatFDColoring color = (MatFDColoring) ctx;
36: MatFDColoringApplyTS(*B,color,t,x1,flag,ts);
37:
38: if (*J != *B) {
39: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
40: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
41: }
42: return(0);
43: }
47: /*@C
48: TSDefaultComputeJacobian - Computes the Jacobian using finite differences.
50: Input Parameters:
51: + ts - TS context
52: . xx1 - compute Jacobian at this point
53: - ctx - application's function context, as set with SNESSetFunction()
55: Output Parameters:
56: + J - Jacobian
57: . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
58: - flag - matrix flag
60: Notes:
61: This routine is slow and expensive, and is not optimized.
63: Sparse approximations using colorings are also available and
64: would be a much better alternative!
66: Level: intermediate
68: .seealso: TSDefaultComputeJacobianColor()
69: @*/
70: PetscErrorCode TSDefaultComputeJacobian(TS ts,PetscReal t,Vec xx1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
71: {
72: Vec f1,f2,xx2;
74: PetscInt i,N,start,end,j;
75: PetscScalar dx,*y,*xx,wscale;
76: PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
77: PetscReal dx_min = 1.e-16,dx_par = 1.e-1;
78: MPI_Comm comm;
79: PetscBool assembled;
80: PetscMPIInt size;
81: const PetscInt *ranges;
82: PetscMPIInt root;
85: VecDuplicate(xx1,&f1);
86: VecDuplicate(xx1,&f2);
87: VecDuplicate(xx1,&xx2);
89: PetscObjectGetComm((PetscObject)xx1,&comm);
90: MPI_Comm_size(comm,&size);
91: MatAssembled(*B,&assembled);
92: if (assembled) {
93: MatZeroEntries(*B);
94: }
96: VecGetSize(xx1,&N);
97: VecGetOwnershipRange(xx1,&start,&end);
98: TSComputeRHSFunction(ts,ts->ptime,xx1,f1);
100: /* Compute Jacobian approximation, 1 column at a time.
101: xx1 = current iterate, f1 = F(xx1)
102: xx2 = perturbed iterate, f2 = F(xx2)
103: */
104: for (i=0; i<N; i++) {
105: VecCopy(xx1,xx2);
106: if (i>= start && i<end) {
107: VecGetArray(xx1,&xx);
108: dx = xx[i-start];
109: VecRestoreArray(xx1,&xx);
110: #if !defined(PETSC_USE_COMPLEX)
111: if (dx < dx_min && dx >= 0.0) dx = dx_par;
112: else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
113: #else
114: if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
115: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
116: #endif
117: dx *= epsilon;
118: wscale = 1.0/dx;
119: VecSetValues(xx2,1,&i,&dx,ADD_VALUES);
120: } else {
121: wscale = 0.0;
122: }
123: TSComputeRHSFunction(ts,t,xx2,f2);
124: VecAXPY(f2,-1.0,f1);
125: /* Communicate scale=1/dx_i to all processors */
126: VecGetOwnershipRanges(xx1,&ranges);
127: root = size;
128: for (j=size-1; j>-1; j--){
129: root--;
130: if (i>=ranges[j]) break;
131: }
132: MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);
134: VecScale(f2,wscale);
135: VecNorm(f2,NORM_INFINITY,&amax); amax *= 1.e-14;
136: VecGetArray(f2,&y);
137: for (j=start; j<end; j++) {
138: if (PetscAbsScalar(y[j-start]) > amax || j == i) {
139: MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);
140: }
141: }
142: VecRestoreArray(f2,&y);
143: }
144: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
145: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
146: if (*B != *J) {
147: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
148: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
149: }
150: *flag = DIFFERENT_NONZERO_PATTERN;
152: VecDestroy(&f1);
153: VecDestroy(&f2);
154: VecDestroy(&xx2);
155: return(0);
156: }