xref: /petsc/src/snes/impls/vi/ss/viss.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3)
1c2fc9fa9SBarry Smith 
2c2fc9fa9SBarry Smith #include <../src/snes/impls/vi/ss/vissimpl.h> /*I "petscsnes.h" I*/
3c2fc9fa9SBarry Smith 
4c2fc9fa9SBarry Smith /*
5c2fc9fa9SBarry Smith   SNESVIComputeMeritFunction - Evaluates the merit function for the mixed complementarity problem.
6c2fc9fa9SBarry Smith 
7c2fc9fa9SBarry Smith   Input Parameter:
8c2fc9fa9SBarry Smith . phi - the semismooth function
9c2fc9fa9SBarry Smith 
10c2fc9fa9SBarry Smith   Output Parameter:
11c2fc9fa9SBarry Smith . merit - the merit function
12c2fc9fa9SBarry Smith . phinorm - ||phi||
13c2fc9fa9SBarry Smith 
14c2fc9fa9SBarry Smith   Notes:
15c2fc9fa9SBarry Smith   The merit function for the mixed complementarity problem is defined as
16c2fc9fa9SBarry Smith      merit = 0.5*phi^T*phi
17c2fc9fa9SBarry Smith */
18c2fc9fa9SBarry Smith #undef __FUNCT__
19c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeMeritFunction"
20c2fc9fa9SBarry Smith static PetscErrorCode SNESVIComputeMeritFunction(Vec phi, PetscReal *merit,PetscReal *phinorm)
21c2fc9fa9SBarry Smith {
22c2fc9fa9SBarry Smith   PetscErrorCode ierr;
23c2fc9fa9SBarry Smith 
24c2fc9fa9SBarry Smith   PetscFunctionBegin;
25c2fc9fa9SBarry Smith   ierr = VecNormBegin(phi,NORM_2,phinorm);CHKERRQ(ierr);
26c2fc9fa9SBarry Smith   ierr = VecNormEnd(phi,NORM_2,phinorm);CHKERRQ(ierr);
27c2fc9fa9SBarry Smith 
28c2fc9fa9SBarry Smith   *merit = 0.5*(*phinorm)*(*phinorm);
29c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
30c2fc9fa9SBarry Smith }
31c2fc9fa9SBarry Smith 
32c2fc9fa9SBarry Smith PETSC_STATIC_INLINE PetscScalar Phi(PetscScalar a,PetscScalar b)
33c2fc9fa9SBarry Smith {
34c2fc9fa9SBarry Smith   return a + b - PetscSqrtScalar(a*a + b*b);
35c2fc9fa9SBarry Smith }
36c2fc9fa9SBarry Smith 
37c2fc9fa9SBarry Smith PETSC_STATIC_INLINE PetscScalar DPhi(PetscScalar a,PetscScalar b)
38c2fc9fa9SBarry Smith {
39c2fc9fa9SBarry Smith   if ((PetscAbsScalar(a) >= 1.e-6) || (PetscAbsScalar(b) >= 1.e-6)) return 1.0 - a/ PetscSqrtScalar(a*a + b*b);
40c2fc9fa9SBarry Smith   else return .5;
41c2fc9fa9SBarry Smith }
42c2fc9fa9SBarry Smith 
43c2fc9fa9SBarry Smith /*
44c2fc9fa9SBarry Smith    SNESVIComputeFunction - Reformulates a system of nonlinear equations in mixed complementarity form to a system of nonlinear equations in semismooth form.
45c2fc9fa9SBarry Smith 
46c2fc9fa9SBarry Smith    Input Parameters:
47c2fc9fa9SBarry Smith .  snes - the SNES context
48d5f1b7e6SEd Bueler .  X - current iterate
49c2fc9fa9SBarry Smith .  functx - user defined function context
50c2fc9fa9SBarry Smith 
51c2fc9fa9SBarry Smith    Output Parameters:
52c2fc9fa9SBarry Smith .  phi - Semismooth function
53c2fc9fa9SBarry Smith 
54c2fc9fa9SBarry Smith */
55c2fc9fa9SBarry Smith #undef __FUNCT__
56c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeFunction"
57c2fc9fa9SBarry Smith static PetscErrorCode SNESVIComputeFunction(SNES snes,Vec X,Vec phi,void *functx)
58c2fc9fa9SBarry Smith {
59c2fc9fa9SBarry Smith   PetscErrorCode    ierr;
60f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*)snes->data;
61c2fc9fa9SBarry Smith   Vec               Xl  = snes->xl,Xu = snes->xu,F = snes->vec_func;
625edff71fSBarry Smith   PetscScalar       *phi_arr,*f_arr,*l,*u;
635edff71fSBarry Smith   const PetscScalar *x_arr;
64c2fc9fa9SBarry Smith   PetscInt          i,nlocal;
65c2fc9fa9SBarry Smith 
66c2fc9fa9SBarry Smith   PetscFunctionBegin;
67c2fc9fa9SBarry Smith   ierr = (*vi->computeuserfunction)(snes,X,F,functx);CHKERRQ(ierr);
68c2fc9fa9SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
695edff71fSBarry Smith   ierr = VecGetArrayRead(X,&x_arr);CHKERRQ(ierr);
70c2fc9fa9SBarry Smith   ierr = VecGetArray(F,&f_arr);CHKERRQ(ierr);
71c2fc9fa9SBarry Smith   ierr = VecGetArray(Xl,&l);CHKERRQ(ierr);
72c2fc9fa9SBarry Smith   ierr = VecGetArray(Xu,&u);CHKERRQ(ierr);
73c2fc9fa9SBarry Smith   ierr = VecGetArray(phi,&phi_arr);CHKERRQ(ierr);
74c2fc9fa9SBarry Smith 
75c2fc9fa9SBarry Smith   for (i=0; i < nlocal; i++) {
76e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
77c2fc9fa9SBarry Smith       phi_arr[i] = f_arr[i];
78e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {                      /* upper bound on variable only */
79c2fc9fa9SBarry Smith       phi_arr[i] = -Phi(u[i] - x_arr[i],-f_arr[i]);
80e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {                       /* lower bound on variable only */
81c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],f_arr[i]);
82c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) {
83c2fc9fa9SBarry Smith       phi_arr[i] = l[i] - x_arr[i];
84c2fc9fa9SBarry Smith     } else {                                                /* both bounds on variable */
85c2fc9fa9SBarry Smith       phi_arr[i] = Phi(x_arr[i] - l[i],-Phi(u[i] - x_arr[i],-f_arr[i]));
86c2fc9fa9SBarry Smith     }
87c2fc9fa9SBarry Smith   }
88c2fc9fa9SBarry Smith 
895edff71fSBarry Smith   ierr = VecRestoreArrayRead(X,&x_arr);CHKERRQ(ierr);
90c2fc9fa9SBarry Smith   ierr = VecRestoreArray(F,&f_arr);CHKERRQ(ierr);
91c2fc9fa9SBarry Smith   ierr = VecRestoreArray(Xl,&l);CHKERRQ(ierr);
92c2fc9fa9SBarry Smith   ierr = VecRestoreArray(Xu,&u);CHKERRQ(ierr);
93c2fc9fa9SBarry Smith   ierr = VecRestoreArray(phi,&phi_arr);CHKERRQ(ierr);
94c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
95c2fc9fa9SBarry Smith }
96c2fc9fa9SBarry Smith 
97c2fc9fa9SBarry Smith /*
98c2fc9fa9SBarry Smith    SNESVIComputeBsubdifferentialVectors - Computes the diagonal shift (Da) and row scaling (Db) vectors needed for the
99c2fc9fa9SBarry Smith                                           the semismooth jacobian.
100c2fc9fa9SBarry Smith */
101c2fc9fa9SBarry Smith #undef __FUNCT__
102c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeBsubdifferentialVectors"
103c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeBsubdifferentialVectors(SNES snes,Vec X,Vec F,Mat jac,Vec Da,Vec Db)
104c2fc9fa9SBarry Smith {
105c2fc9fa9SBarry Smith   PetscErrorCode ierr;
106c2fc9fa9SBarry Smith   PetscScalar    *l,*u,*x,*f,*da,*db,da1,da2,db1,db2;
107c2fc9fa9SBarry Smith   PetscInt       i,nlocal;
108c2fc9fa9SBarry Smith 
109c2fc9fa9SBarry Smith   PetscFunctionBegin;
110c2fc9fa9SBarry Smith   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
111c2fc9fa9SBarry Smith   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
112c2fc9fa9SBarry Smith   ierr = VecGetArray(snes->xl,&l);CHKERRQ(ierr);
113c2fc9fa9SBarry Smith   ierr = VecGetArray(snes->xu,&u);CHKERRQ(ierr);
114c2fc9fa9SBarry Smith   ierr = VecGetArray(Da,&da);CHKERRQ(ierr);
115c2fc9fa9SBarry Smith   ierr = VecGetArray(Db,&db);CHKERRQ(ierr);
116c2fc9fa9SBarry Smith   ierr = VecGetLocalSize(X,&nlocal);CHKERRQ(ierr);
117c2fc9fa9SBarry Smith 
118c2fc9fa9SBarry Smith   for (i=0; i< nlocal; i++) {
119e270355aSBarry Smith     if ((PetscRealPart(l[i]) <= PETSC_NINFINITY) && (PetscRealPart(u[i]) >= PETSC_INFINITY)) { /* no constraints on variable */
120c2fc9fa9SBarry Smith       da[i] = 0;
121c2fc9fa9SBarry Smith       db[i] = 1;
122e270355aSBarry Smith     } else if (PetscRealPart(l[i]) <= PETSC_NINFINITY) {                     /* upper bound on variable only */
123c2fc9fa9SBarry Smith       da[i] = DPhi(u[i] - x[i], -f[i]);
124c2fc9fa9SBarry Smith       db[i] = DPhi(-f[i],u[i] - x[i]);
125e270355aSBarry Smith     } else if (PetscRealPart(u[i]) >= PETSC_INFINITY) {                      /* lower bound on variable only */
126c2fc9fa9SBarry Smith       da[i] = DPhi(x[i] - l[i], f[i]);
127c2fc9fa9SBarry Smith       db[i] = DPhi(f[i],x[i] - l[i]);
128c2fc9fa9SBarry Smith     } else if (l[i] == u[i]) {                              /* fixed variable */
129c2fc9fa9SBarry Smith       da[i] = 1;
130c2fc9fa9SBarry Smith       db[i] = 0;
131c2fc9fa9SBarry Smith     } else {                                                /* upper and lower bounds on variable */
132c2fc9fa9SBarry Smith       da1   = DPhi(x[i] - l[i], -Phi(u[i] - x[i], -f[i]));
133c2fc9fa9SBarry Smith       db1   = DPhi(-Phi(u[i] - x[i], -f[i]),x[i] - l[i]);
134c2fc9fa9SBarry Smith       da2   = DPhi(u[i] - x[i], -f[i]);
135c2fc9fa9SBarry Smith       db2   = DPhi(-f[i],u[i] - x[i]);
136c2fc9fa9SBarry Smith       da[i] = da1 + db1*da2;
137c2fc9fa9SBarry Smith       db[i] = db1*db2;
138c2fc9fa9SBarry Smith     }
139c2fc9fa9SBarry Smith   }
140c2fc9fa9SBarry Smith 
141c2fc9fa9SBarry Smith   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
142c2fc9fa9SBarry Smith   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
143c2fc9fa9SBarry Smith   ierr = VecRestoreArray(snes->xl,&l);CHKERRQ(ierr);
144c2fc9fa9SBarry Smith   ierr = VecRestoreArray(snes->xu,&u);CHKERRQ(ierr);
145c2fc9fa9SBarry Smith   ierr = VecRestoreArray(Da,&da);CHKERRQ(ierr);
146c2fc9fa9SBarry Smith   ierr = VecRestoreArray(Db,&db);CHKERRQ(ierr);
147c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
148c2fc9fa9SBarry Smith }
149c2fc9fa9SBarry Smith 
150c2fc9fa9SBarry Smith /*
151c2fc9fa9SBarry Smith    SNESVIComputeJacobian - Computes the jacobian of the semismooth function.The Jacobian for the semismooth function is an element of the B-subdifferential of the Fischer-Burmeister function for complementarity problems.
152c2fc9fa9SBarry Smith 
153c2fc9fa9SBarry Smith    Input Parameters:
154c2fc9fa9SBarry Smith .  Da       - Diagonal shift vector for the semismooth jacobian.
155c2fc9fa9SBarry Smith .  Db       - Row scaling vector for the semismooth jacobian.
156c2fc9fa9SBarry Smith 
157c2fc9fa9SBarry Smith    Output Parameters:
158c2fc9fa9SBarry Smith .  jac      - semismooth jacobian
159c2fc9fa9SBarry Smith .  jac_pre  - optional preconditioning matrix
160c2fc9fa9SBarry Smith 
161c2fc9fa9SBarry Smith    Notes:
162c2fc9fa9SBarry Smith    The semismooth jacobian matrix is given by
163c2fc9fa9SBarry Smith    jac = Da + Db*jacfun
164c2fc9fa9SBarry Smith    where Db is the row scaling matrix stored as a vector,
165c2fc9fa9SBarry Smith          Da is the diagonal perturbation matrix stored as a vector
166c2fc9fa9SBarry Smith    and   jacfun is the jacobian of the original nonlinear function.
167c2fc9fa9SBarry Smith */
168c2fc9fa9SBarry Smith #undef __FUNCT__
169c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeJacobian"
170c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeJacobian(Mat jac, Mat jac_pre,Vec Da, Vec Db)
171c2fc9fa9SBarry Smith {
172c2fc9fa9SBarry Smith   PetscErrorCode ierr;
173c2fc9fa9SBarry Smith 
174c2fc9fa9SBarry Smith   /* Do row scaling  and add diagonal perturbation */
1750298fd71SBarry Smith   ierr = MatDiagonalScale(jac,Db,NULL);CHKERRQ(ierr);
176c2fc9fa9SBarry Smith   ierr = MatDiagonalSet(jac,Da,ADD_VALUES);CHKERRQ(ierr);
177c2fc9fa9SBarry Smith   if (jac != jac_pre) { /* If jac and jac_pre are different */
1780298fd71SBarry Smith     ierr = MatDiagonalScale(jac_pre,Db,NULL);CHKERRQ(ierr);
179c2fc9fa9SBarry Smith     ierr = MatDiagonalSet(jac_pre,Da,ADD_VALUES);CHKERRQ(ierr);
180c2fc9fa9SBarry Smith   }
181c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
182c2fc9fa9SBarry Smith }
183c2fc9fa9SBarry Smith 
184c2fc9fa9SBarry Smith /*
185c2fc9fa9SBarry Smith    SNESVIComputeMeritFunctionGradient - Computes the gradient of the merit function psi.
186c2fc9fa9SBarry Smith 
187c2fc9fa9SBarry Smith    Input Parameters:
188c2fc9fa9SBarry Smith    phi - semismooth function.
189c2fc9fa9SBarry Smith    H   - semismooth jacobian
190c2fc9fa9SBarry Smith 
191c2fc9fa9SBarry Smith    Output Parameters:
192c2fc9fa9SBarry Smith    dpsi - merit function gradient
193c2fc9fa9SBarry Smith 
194c2fc9fa9SBarry Smith    Notes:
195c2fc9fa9SBarry Smith   The merit function gradient is computed as follows
196c2fc9fa9SBarry Smith         dpsi = H^T*phi
197c2fc9fa9SBarry Smith */
198c2fc9fa9SBarry Smith #undef __FUNCT__
199c2fc9fa9SBarry Smith #define __FUNCT__ "SNESVIComputeMeritFunctionGradient"
200c2fc9fa9SBarry Smith PetscErrorCode SNESVIComputeMeritFunctionGradient(Mat H, Vec phi, Vec dpsi)
201c2fc9fa9SBarry Smith {
202c2fc9fa9SBarry Smith   PetscErrorCode ierr;
203c2fc9fa9SBarry Smith 
204c2fc9fa9SBarry Smith   PetscFunctionBegin;
205c2fc9fa9SBarry Smith   ierr = MatMultTranspose(H,phi,dpsi);CHKERRQ(ierr);
206c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
207c2fc9fa9SBarry Smith }
208c2fc9fa9SBarry Smith 
209c2fc9fa9SBarry Smith 
210c2fc9fa9SBarry Smith 
211c2fc9fa9SBarry Smith /*
212f450aa47SBarry Smith    SNESSolve_VINEWTONSSLS - Solves the complementarity problem with a semismooth Newton
213c2fc9fa9SBarry Smith    method using a line search.
214c2fc9fa9SBarry Smith 
215c2fc9fa9SBarry Smith    Input Parameters:
216c2fc9fa9SBarry Smith .  snes - the SNES context
217c2fc9fa9SBarry Smith 
218c2fc9fa9SBarry Smith    Application Interface Routine: SNESSolve()
219c2fc9fa9SBarry Smith 
220c2fc9fa9SBarry Smith    Notes:
221c2fc9fa9SBarry Smith    This implements essentially a semismooth Newton method with a
222d5f1b7e6SEd Bueler    line search. The default line search does not do any line search
223d5f1b7e6SEd Bueler    but rather takes a full Newton step.
224016d19f9SBarry Smith 
22504d7464bSBarry Smith    Developer Note: the code in this file should be slightly modified so that this routine need not exist and the SNESSolve_NEWTONLS() routine is called directly with the appropriate wrapped function and Jacobian evaluations
226016d19f9SBarry Smith 
227c2fc9fa9SBarry Smith */
228c2fc9fa9SBarry Smith #undef __FUNCT__
229f450aa47SBarry Smith #define __FUNCT__ "SNESSolve_VINEWTONSSLS"
230f450aa47SBarry Smith PetscErrorCode SNESSolve_VINEWTONSSLS(SNES snes)
231c2fc9fa9SBarry Smith {
232f450aa47SBarry Smith   SNES_VINEWTONSSLS    *vi = (SNES_VINEWTONSSLS*)snes->data;
233c2fc9fa9SBarry Smith   PetscErrorCode       ierr;
234c2fc9fa9SBarry Smith   PetscInt             maxits,i,lits;
235422a814eSBarry Smith   SNESLineSearchReason lssucceed;
236c2fc9fa9SBarry Smith   PetscReal            gnorm,xnorm=0,ynorm;
2379bd66eb0SPeter Brune   Vec                  Y,X,F;
238c2fc9fa9SBarry Smith   KSPConvergedReason   kspreason;
2396cab3a1bSJed Brown   DM                   dm;
240942e3340SBarry Smith   DMSNES               sdm;
241c2fc9fa9SBarry Smith 
242c2fc9fa9SBarry Smith   PetscFunctionBegin;
2436cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
244942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
2451aa26658SKarl Rupp 
24622c6f798SBarry Smith   vi->computeuserfunction   = sdm->ops->computefunction;
24722c6f798SBarry Smith   sdm->ops->computefunction = SNESVIComputeFunction;
248c2fc9fa9SBarry Smith 
249c2fc9fa9SBarry Smith   snes->numFailures            = 0;
250c2fc9fa9SBarry Smith   snes->numLinearSolveFailures = 0;
251c2fc9fa9SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
252c2fc9fa9SBarry Smith 
253c2fc9fa9SBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
254c2fc9fa9SBarry Smith   X      = snes->vec_sol;               /* solution vector */
255c2fc9fa9SBarry Smith   F      = snes->vec_func;              /* residual vector */
256c2fc9fa9SBarry Smith   Y      = snes->work[0];               /* work vectors */
257c2fc9fa9SBarry Smith 
258e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
259c2fc9fa9SBarry Smith   snes->iter = 0;
260c2fc9fa9SBarry Smith   snes->norm = 0.0;
261e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
262c2fc9fa9SBarry Smith 
263c2fc9fa9SBarry Smith   ierr = SNESVIProjectOntoBounds(snes,X);CHKERRQ(ierr);
264c2fc9fa9SBarry Smith   ierr = SNESComputeFunction(snes,X,vi->phi);CHKERRQ(ierr);
265c2fc9fa9SBarry Smith   if (snes->domainerror) {
266c2fc9fa9SBarry Smith     snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
26722c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
268c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
269c2fc9fa9SBarry Smith   }
270c2fc9fa9SBarry Smith   /* Compute Merit function */
271c2fc9fa9SBarry Smith   ierr = SNESVIComputeMeritFunction(vi->phi,&vi->merit,&vi->phinorm);CHKERRQ(ierr);
272c2fc9fa9SBarry Smith 
273c2fc9fa9SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);        /* xnorm <- ||x||  */
274c2fc9fa9SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
275422a814eSBarry Smith   SNESCheckFunctionNorm(snes,vi->merit);
276c2fc9fa9SBarry Smith 
277e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
278c2fc9fa9SBarry Smith   snes->norm = vi->phinorm;
279e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
280a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,vi->phinorm,0);CHKERRQ(ierr);
281c2fc9fa9SBarry Smith   ierr       = SNESMonitor(snes,0,vi->phinorm);CHKERRQ(ierr);
282c2fc9fa9SBarry Smith 
283c2fc9fa9SBarry Smith   /* test convergence */
284c2fc9fa9SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
285c2fc9fa9SBarry Smith   if (snes->reason) {
28622c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
287c2fc9fa9SBarry Smith     PetscFunctionReturn(0);
288c2fc9fa9SBarry Smith   }
289c2fc9fa9SBarry Smith 
290c2fc9fa9SBarry Smith   for (i=0; i<maxits; i++) {
291c2fc9fa9SBarry Smith 
292c2fc9fa9SBarry Smith     /* Call general purpose update function */
293c2fc9fa9SBarry Smith     if (snes->ops->update) {
294c2fc9fa9SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
295c2fc9fa9SBarry Smith     }
296c2fc9fa9SBarry Smith 
297c2fc9fa9SBarry Smith     /* Solve J Y = Phi, where J is the semismooth jacobian */
298b9600128SPeter Brune 
299b9600128SPeter Brune     /* Get the jacobian -- note that the function must be the original function for snes_fd and snes_fd_color to work for this*/
30022c6f798SBarry Smith     sdm->ops->computefunction = vi->computeuserfunction;
301d1e9a80fSBarry Smith     ierr                      = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
30222c6f798SBarry Smith     sdm->ops->computefunction = SNESVIComputeFunction;
303b9600128SPeter Brune 
304c2fc9fa9SBarry Smith     /* Get the diagonal shift and row scaling vectors */
305c2fc9fa9SBarry Smith     ierr = SNESVIComputeBsubdifferentialVectors(snes,X,F,snes->jacobian,vi->Da,vi->Db);CHKERRQ(ierr);
306c2fc9fa9SBarry Smith     /* Compute the semismooth jacobian */
307c2fc9fa9SBarry Smith     ierr = SNESVIComputeJacobian(snes->jacobian,snes->jacobian_pre,vi->Da,vi->Db);CHKERRQ(ierr);
308c2fc9fa9SBarry Smith     /* Compute the merit function gradient */
309c2fc9fa9SBarry Smith     ierr = SNESVIComputeMeritFunctionGradient(snes->jacobian,vi->phi,vi->dpsi);CHKERRQ(ierr);
31023ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
311d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,vi->phi,Y);CHKERRQ(ierr);
312c2fc9fa9SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
313c2fc9fa9SBarry Smith 
314c2fc9fa9SBarry Smith     if (kspreason < 0) {
315c2fc9fa9SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
316c2fc9fa9SBarry Smith         ierr         = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
317c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
318c2fc9fa9SBarry Smith         break;
319c2fc9fa9SBarry Smith       }
320c2fc9fa9SBarry Smith     }
321c2fc9fa9SBarry Smith     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
322c2fc9fa9SBarry Smith     snes->linear_its += lits;
323c2fc9fa9SBarry Smith     ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
324c2fc9fa9SBarry Smith     /*
3256b2b7091SBarry Smith     if (snes->ops->precheck) {
326c2fc9fa9SBarry Smith       PetscBool changed_y = PETSC_FALSE;
3276b2b7091SBarry Smith       ierr = (*snes->ops->precheck)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
328c2fc9fa9SBarry Smith     }
329c2fc9fa9SBarry Smith 
330c2fc9fa9SBarry Smith     if (PetscLogPrintInfo) {
331c2fc9fa9SBarry Smith       ierr = SNESVICheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
332c2fc9fa9SBarry Smith     }
333c2fc9fa9SBarry Smith     */
334c2fc9fa9SBarry Smith     /* Compute a (scaled) negative update in the line search routine:
335c2fc9fa9SBarry Smith          Y <- X - lambda*Y
336c2fc9fa9SBarry Smith        and evaluate G = function(Y) (depends on the line search).
337c2fc9fa9SBarry Smith     */
338c2fc9fa9SBarry Smith     ierr  = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
339c2fc9fa9SBarry Smith     ynorm = 1; gnorm = vi->phinorm;
340f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, vi->phi, &gnorm, Y);CHKERRQ(ierr);
341422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr);
342f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
343c2fc9fa9SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)vi->phinorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
344c2fc9fa9SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
345c2fc9fa9SBarry Smith     if (snes->domainerror) {
346c2fc9fa9SBarry Smith       snes->reason              = SNES_DIVERGED_FUNCTION_DOMAIN;
34722c6f798SBarry Smith       sdm->ops->computefunction = vi->computeuserfunction;
348c2fc9fa9SBarry Smith       PetscFunctionReturn(0);
349c2fc9fa9SBarry Smith     }
350422a814eSBarry Smith     if (lssucceed) {
351c2fc9fa9SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
352c2fc9fa9SBarry Smith         PetscBool ismin;
353c2fc9fa9SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3549bd66eb0SPeter Brune         ierr         = SNESVICheckLocalMin_Private(snes,snes->jacobian,vi->phi,X,gnorm,&ismin);CHKERRQ(ierr);
355c2fc9fa9SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
356c2fc9fa9SBarry Smith         break;
357c2fc9fa9SBarry Smith       }
358c2fc9fa9SBarry Smith     }
359c2fc9fa9SBarry Smith     /* Update function and solution vectors */
360c2fc9fa9SBarry Smith     vi->phinorm = gnorm;
361c2fc9fa9SBarry Smith     vi->merit   = 0.5*vi->phinorm*vi->phinorm;
362c2fc9fa9SBarry Smith     /* Monitor convergence */
363e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
364c2fc9fa9SBarry Smith     snes->iter = i+1;
365c2fc9fa9SBarry Smith     snes->norm = vi->phinorm;
366e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
367a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
368c2fc9fa9SBarry Smith     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
369c2fc9fa9SBarry Smith     /* Test for convergence, xnorm = || X || */
370e2a6519dSDmitry Karpeev     if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
371c2fc9fa9SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,vi->phinorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
372c2fc9fa9SBarry Smith     if (snes->reason) break;
373c2fc9fa9SBarry Smith   }
374c2fc9fa9SBarry Smith   if (i == maxits) {
375c2fc9fa9SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
376c2fc9fa9SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
377c2fc9fa9SBarry Smith   }
37822c6f798SBarry Smith   sdm->ops->computefunction = vi->computeuserfunction;
379c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
380c2fc9fa9SBarry Smith }
381c2fc9fa9SBarry Smith 
382c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
383c2fc9fa9SBarry Smith /*
384f450aa47SBarry Smith    SNESSetUp_VINEWTONSSLS - Sets up the internal data structures for the later use
385c2fc9fa9SBarry Smith    of the SNES nonlinear solver.
386c2fc9fa9SBarry Smith 
387c2fc9fa9SBarry Smith    Input Parameter:
388c2fc9fa9SBarry Smith .  snes - the SNES context
389c2fc9fa9SBarry Smith 
390c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetUp()
391c2fc9fa9SBarry Smith 
392c2fc9fa9SBarry Smith    Notes:
393c2fc9fa9SBarry Smith    For basic use of the SNES solvers, the user need not explicitly call
394c2fc9fa9SBarry Smith    SNESSetUp(), since these actions will automatically occur during
395c2fc9fa9SBarry Smith    the call to SNESSolve().
396c2fc9fa9SBarry Smith  */
397c2fc9fa9SBarry Smith #undef __FUNCT__
398f450aa47SBarry Smith #define __FUNCT__ "SNESSetUp_VINEWTONSSLS"
399f450aa47SBarry Smith PetscErrorCode SNESSetUp_VINEWTONSSLS(SNES snes)
400c2fc9fa9SBarry Smith {
401c2fc9fa9SBarry Smith   PetscErrorCode    ierr;
402f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data;
403c2fc9fa9SBarry Smith 
404c2fc9fa9SBarry Smith   PetscFunctionBegin;
405c2fc9fa9SBarry Smith   ierr = SNESSetUp_VI(snes);CHKERRQ(ierr);
406c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->dpsi);CHKERRQ(ierr);
407c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->phi);CHKERRQ(ierr);
408c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->Da);CHKERRQ(ierr);
409c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->Db);CHKERRQ(ierr);
410c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->z);CHKERRQ(ierr);
411c2fc9fa9SBarry Smith   ierr = VecDuplicate(snes->vec_sol, &vi->t);CHKERRQ(ierr);
412c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
413c2fc9fa9SBarry Smith }
414c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
415c2fc9fa9SBarry Smith #undef __FUNCT__
416f450aa47SBarry Smith #define __FUNCT__ "SNESReset_VINEWTONSSLS"
417f450aa47SBarry Smith PetscErrorCode SNESReset_VINEWTONSSLS(SNES snes)
418c2fc9fa9SBarry Smith {
419f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi = (SNES_VINEWTONSSLS*) snes->data;
420c2fc9fa9SBarry Smith   PetscErrorCode    ierr;
421c2fc9fa9SBarry Smith 
422c2fc9fa9SBarry Smith   PetscFunctionBegin;
423c2fc9fa9SBarry Smith   ierr = SNESReset_VI(snes);CHKERRQ(ierr);
424c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->dpsi);CHKERRQ(ierr);
425c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->phi);CHKERRQ(ierr);
426c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->Da);CHKERRQ(ierr);
427c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->Db);CHKERRQ(ierr);
428c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->z);CHKERRQ(ierr);
429c2fc9fa9SBarry Smith   ierr = VecDestroy(&vi->t);CHKERRQ(ierr);
430c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
431c2fc9fa9SBarry Smith }
432c2fc9fa9SBarry Smith 
433c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
434c2fc9fa9SBarry Smith /*
435f450aa47SBarry Smith    SNESSetFromOptions_VINEWTONSSLS - Sets various parameters for the SNESVI method.
436c2fc9fa9SBarry Smith 
437c2fc9fa9SBarry Smith    Input Parameter:
438c2fc9fa9SBarry Smith .  snes - the SNES context
439c2fc9fa9SBarry Smith 
440c2fc9fa9SBarry Smith    Application Interface Routine: SNESSetFromOptions()
441c2fc9fa9SBarry Smith */
442c2fc9fa9SBarry Smith #undef __FUNCT__
443f450aa47SBarry Smith #define __FUNCT__ "SNESSetFromOptions_VINEWTONSSLS"
4444416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_VINEWTONSSLS(PetscOptionItems *PetscOptionsObject,SNES snes)
445c2fc9fa9SBarry Smith {
446c2fc9fa9SBarry Smith   PetscErrorCode ierr;
447f1c6b773SPeter Brune   SNESLineSearch linesearch;
448c2fc9fa9SBarry Smith 
449c2fc9fa9SBarry Smith   PetscFunctionBegin;
450e55864a3SBarry Smith   ierr = SNESSetFromOptions_VI(PetscOptionsObject,snes);CHKERRQ(ierr);
451e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES semismooth method options");CHKERRQ(ierr);
452c2fc9fa9SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4539bd66eb0SPeter Brune   /* set up the default line search */
4549bd66eb0SPeter Brune   if (!snes->linesearch) {
4557601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
4561a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
457639ea3faSPeter Brune     ierr = SNESLineSearchBTSetAlpha(linesearch, 0.0);CHKERRQ(ierr);
4589bd66eb0SPeter Brune   }
459c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
460c2fc9fa9SBarry Smith }
461c2fc9fa9SBarry Smith 
462c2fc9fa9SBarry Smith 
463c2fc9fa9SBarry Smith /* -------------------------------------------------------------------------- */
464c2fc9fa9SBarry Smith /*MC
465f450aa47SBarry Smith       SNESVINEWTONSSLS - Semi-smooth solver for variational inequalities based on Newton's method
466c2fc9fa9SBarry Smith 
46761589011SJed Brown    Options Database:
468d5f1b7e6SEd Bueler +   -snes_type <vinewtonssls,vinewtonrsls> a semi-smooth solver, a reduced space active set method
46961589011SJed Brown -   -snes_vi_monitor - prints the number of active constraints at each iteration.
47061589011SJed Brown 
471c2fc9fa9SBarry Smith    Level: beginner
472c2fc9fa9SBarry Smith 
473b80f3ac1SShri Abhyankar    References:
47496a0c994SBarry Smith +  1. -  T. S. Munson, F. Facchinei, M. C. Ferris, A. Fischer, and C. Kanzow. The semismooth
475b80f3ac1SShri Abhyankar      algorithm for large scale complementarity problems. INFORMS Journal on Computing, 13 (2001).
47696a0c994SBarry Smith -  2. -  T. S. Munson, and S. Benson. Flexible Complementarity Solvers for Large Scale
477d5f1b7e6SEd Bueler      Applications, Optimization Methods and Software, 21 (2006).
478b80f3ac1SShri Abhyankar 
479d5f1b7e6SEd Bueler .seealso:  SNESVISetVariableBounds(), SNESVISetComputeVariableBounds(), SNESCreate(), SNES, SNESSetType(), SNESVINEWTONRSLS, SNESNEWTONTR, SNESLineSearchSet(),SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
480c2fc9fa9SBarry Smith 
481c2fc9fa9SBarry Smith M*/
482c2fc9fa9SBarry Smith #undef __FUNCT__
483f450aa47SBarry Smith #define __FUNCT__ "SNESCreate_VINEWTONSSLS"
4848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES snes)
485c2fc9fa9SBarry Smith {
486c2fc9fa9SBarry Smith   PetscErrorCode    ierr;
487f450aa47SBarry Smith   SNES_VINEWTONSSLS *vi;
488c2fc9fa9SBarry Smith 
489c2fc9fa9SBarry Smith   PetscFunctionBegin;
490f450aa47SBarry Smith   snes->ops->reset          = SNESReset_VINEWTONSSLS;
491f450aa47SBarry Smith   snes->ops->setup          = SNESSetUp_VINEWTONSSLS;
492f450aa47SBarry Smith   snes->ops->solve          = SNESSolve_VINEWTONSSLS;
493c2fc9fa9SBarry Smith   snes->ops->destroy        = SNESDestroy_VI;
494f450aa47SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_VINEWTONSSLS;
4950298fd71SBarry Smith   snes->ops->view           = NULL;
496c2fc9fa9SBarry Smith 
497c2fc9fa9SBarry Smith   snes->usesksp = PETSC_TRUE;
498c2fc9fa9SBarry Smith   snes->usespc  = PETSC_FALSE;
499c2fc9fa9SBarry Smith 
500*4fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
501*4fc747eaSLawrence Mitchell 
502b00a9115SJed Brown   ierr       = PetscNewLog(snes,&vi);CHKERRQ(ierr);
503c2fc9fa9SBarry Smith   snes->data = (void*)vi;
504c2fc9fa9SBarry Smith 
505bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESVISetVariableBounds_C",SNESVISetVariableBounds_VI);CHKERRQ(ierr);
506bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESVISetComputeVariableBounds_C",SNESVISetComputeVariableBounds_VI);CHKERRQ(ierr);
507c2fc9fa9SBarry Smith   PetscFunctionReturn(0);
508c2fc9fa9SBarry Smith }
509c2fc9fa9SBarry Smith 
510