LCOV - code coverage report
Current view: top level - nep/interface - nepbasic.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 371 397 93.5 %
Date: 2024-11-21 00:34:55 Functions: 24 27 88.9 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : /*
      11             :    Basic NEP routines
      12             : */
      13             : 
      14             : #include <slepc/private/nepimpl.h>      /*I "slepcnep.h" I*/
      15             : 
      16             : /* Logging support */
      17             : PetscClassId      NEP_CLASSID = 0;
      18             : PetscLogEvent     NEP_SetUp = 0,NEP_Solve = 0,NEP_Refine = 0,NEP_FunctionEval = 0,NEP_JacobianEval = 0,NEP_Resolvent = 0,NEP_CISS_SVD = 0;
      19             : 
      20             : /* List of registered NEP routines */
      21             : PetscFunctionList NEPList = NULL;
      22             : PetscBool         NEPRegisterAllCalled = PETSC_FALSE;
      23             : 
      24             : /* List of registered NEP monitors */
      25             : PetscFunctionList NEPMonitorList              = NULL;
      26             : PetscFunctionList NEPMonitorCreateList        = NULL;
      27             : PetscFunctionList NEPMonitorDestroyList       = NULL;
      28             : PetscBool         NEPMonitorRegisterAllCalled = PETSC_FALSE;
      29             : 
      30             : /*@
      31             :    NEPCreate - Creates the default NEP context.
      32             : 
      33             :    Collective
      34             : 
      35             :    Input Parameter:
      36             : .  comm - MPI communicator
      37             : 
      38             :    Output Parameter:
      39             : .  outnep - location to put the NEP context
      40             : 
      41             :    Level: beginner
      42             : 
      43             : .seealso: NEPSetUp(), NEPSolve(), NEPDestroy(), NEP
      44             : @*/
      45         102 : PetscErrorCode NEPCreate(MPI_Comm comm,NEP *outnep)
      46             : {
      47         102 :   NEP            nep;
      48             : 
      49         102 :   PetscFunctionBegin;
      50         102 :   PetscAssertPointer(outnep,2);
      51         102 :   PetscCall(NEPInitializePackage());
      52         102 :   PetscCall(SlepcHeaderCreate(nep,NEP_CLASSID,"NEP","Nonlinear Eigenvalue Problem","NEP",comm,NEPDestroy,NEPView));
      53             : 
      54         102 :   nep->max_it          = PETSC_DETERMINE;
      55         102 :   nep->nev             = 1;
      56         102 :   nep->ncv             = PETSC_DETERMINE;
      57         102 :   nep->mpd             = PETSC_DETERMINE;
      58         102 :   nep->nini            = 0;
      59         102 :   nep->target          = 0.0;
      60         102 :   nep->tol             = PETSC_DETERMINE;
      61         102 :   nep->conv            = NEP_CONV_REL;
      62         102 :   nep->stop            = NEP_STOP_BASIC;
      63         102 :   nep->which           = (NEPWhich)0;
      64         102 :   nep->problem_type    = (NEPProblemType)0;
      65         102 :   nep->refine          = NEP_REFINE_NONE;
      66         102 :   nep->npart           = 1;
      67         102 :   nep->rtol            = PETSC_DETERMINE;
      68         102 :   nep->rits            = PETSC_DETERMINE;
      69         102 :   nep->scheme          = (NEPRefineScheme)0;
      70         102 :   nep->trackall        = PETSC_FALSE;
      71         102 :   nep->twosided        = PETSC_FALSE;
      72             : 
      73         102 :   nep->computefunction = NULL;
      74         102 :   nep->computejacobian = NULL;
      75         102 :   nep->functionctx     = NULL;
      76         102 :   nep->jacobianctx     = NULL;
      77         102 :   nep->converged       = NEPConvergedRelative;
      78         102 :   nep->convergeduser   = NULL;
      79         102 :   nep->convergeddestroy= NULL;
      80         102 :   nep->stopping        = NEPStoppingBasic;
      81         102 :   nep->stoppinguser    = NULL;
      82         102 :   nep->stoppingdestroy = NULL;
      83         102 :   nep->convergedctx    = NULL;
      84         102 :   nep->stoppingctx     = NULL;
      85         102 :   nep->numbermonitors  = 0;
      86             : 
      87         102 :   nep->ds              = NULL;
      88         102 :   nep->V               = NULL;
      89         102 :   nep->W               = NULL;
      90         102 :   nep->rg              = NULL;
      91         102 :   nep->function        = NULL;
      92         102 :   nep->function_pre    = NULL;
      93         102 :   nep->jacobian        = NULL;
      94         102 :   nep->A               = NULL;
      95         102 :   nep->f               = NULL;
      96         102 :   nep->nt              = 0;
      97         102 :   nep->mstr            = UNKNOWN_NONZERO_PATTERN;
      98         102 :   nep->P               = NULL;
      99         102 :   nep->mstrp           = UNKNOWN_NONZERO_PATTERN;
     100         102 :   nep->IS              = NULL;
     101         102 :   nep->eigr            = NULL;
     102         102 :   nep->eigi            = NULL;
     103         102 :   nep->errest          = NULL;
     104         102 :   nep->perm            = NULL;
     105         102 :   nep->nwork           = 0;
     106         102 :   nep->work            = NULL;
     107         102 :   nep->data            = NULL;
     108             : 
     109         102 :   nep->state           = NEP_STATE_INITIAL;
     110         102 :   nep->nconv           = 0;
     111         102 :   nep->its             = 0;
     112         102 :   nep->n               = 0;
     113         102 :   nep->nloc            = 0;
     114         102 :   nep->nrma            = NULL;
     115         102 :   nep->fui             = (NEPUserInterface)0;
     116         102 :   nep->useds           = PETSC_FALSE;
     117         102 :   nep->resolvent       = NULL;
     118         102 :   nep->reason          = NEP_CONVERGED_ITERATING;
     119             : 
     120         102 :   PetscCall(PetscNew(&nep->sc));
     121         102 :   *outnep = nep;
     122         102 :   PetscFunctionReturn(PETSC_SUCCESS);
     123             : }
     124             : 
     125             : /*@
     126             :    NEPSetType - Selects the particular solver to be used in the NEP object.
     127             : 
     128             :    Logically Collective
     129             : 
     130             :    Input Parameters:
     131             : +  nep      - the nonlinear eigensolver context
     132             : -  type     - a known method
     133             : 
     134             :    Options Database Key:
     135             : .  -nep_type <method> - Sets the method; use -help for a list
     136             :     of available methods
     137             : 
     138             :    Notes:
     139             :    See "slepc/include/slepcnep.h" for available methods.
     140             : 
     141             :    Normally, it is best to use the NEPSetFromOptions() command and
     142             :    then set the NEP type from the options database rather than by using
     143             :    this routine.  Using the options database provides the user with
     144             :    maximum flexibility in evaluating the different available methods.
     145             :    The NEPSetType() routine is provided for those situations where it
     146             :    is necessary to set the iterative solver independently of the command
     147             :    line or options database.
     148             : 
     149             :    Level: intermediate
     150             : 
     151             : .seealso: NEPType
     152             : @*/
     153         108 : PetscErrorCode NEPSetType(NEP nep,NEPType type)
     154             : {
     155         108 :   PetscErrorCode (*r)(NEP);
     156         108 :   PetscBool      match;
     157             : 
     158         108 :   PetscFunctionBegin;
     159         108 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     160         108 :   PetscAssertPointer(type,2);
     161             : 
     162         108 :   PetscCall(PetscObjectTypeCompare((PetscObject)nep,type,&match));
     163         108 :   if (match) PetscFunctionReturn(PETSC_SUCCESS);
     164             : 
     165         105 :   PetscCall(PetscFunctionListFind(NEPList,type,&r));
     166         105 :   PetscCheck(r,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown NEP type given: %s",type);
     167             : 
     168         105 :   PetscTryTypeMethod(nep,destroy);
     169         105 :   PetscCall(PetscMemzero(nep->ops,sizeof(struct _NEPOps)));
     170             : 
     171         105 :   nep->state = NEP_STATE_INITIAL;
     172         105 :   PetscCall(PetscObjectChangeTypeName((PetscObject)nep,type));
     173         105 :   PetscCall((*r)(nep));
     174         105 :   PetscFunctionReturn(PETSC_SUCCESS);
     175             : }
     176             : 
     177             : /*@
     178             :    NEPGetType - Gets the NEP type as a string from the NEP object.
     179             : 
     180             :    Not Collective
     181             : 
     182             :    Input Parameter:
     183             : .  nep - the eigensolver context
     184             : 
     185             :    Output Parameter:
     186             : .  type - name of NEP method
     187             : 
     188             :    Level: intermediate
     189             : 
     190             : .seealso: NEPSetType()
     191             : @*/
     192          17 : PetscErrorCode NEPGetType(NEP nep,NEPType *type)
     193             : {
     194          17 :   PetscFunctionBegin;
     195          17 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     196          17 :   PetscAssertPointer(type,2);
     197          17 :   *type = ((PetscObject)nep)->type_name;
     198          17 :   PetscFunctionReturn(PETSC_SUCCESS);
     199             : }
     200             : 
     201             : /*@C
     202             :    NEPRegister - Adds a method to the nonlinear eigenproblem solver package.
     203             : 
     204             :    Not Collective
     205             : 
     206             :    Input Parameters:
     207             : +  name - name of a new user-defined solver
     208             : -  function - routine to create the solver context
     209             : 
     210             :    Notes:
     211             :    NEPRegister() may be called multiple times to add several user-defined solvers.
     212             : 
     213             :    Example Usage:
     214             : .vb
     215             :     NEPRegister("my_solver",MySolverCreate);
     216             : .ve
     217             : 
     218             :    Then, your solver can be chosen with the procedural interface via
     219             : $     NEPSetType(nep,"my_solver")
     220             :    or at runtime via the option
     221             : $     -nep_type my_solver
     222             : 
     223             :    Level: advanced
     224             : 
     225             : .seealso: NEPRegisterAll()
     226             : @*/
     227         515 : PetscErrorCode NEPRegister(const char *name,PetscErrorCode (*function)(NEP))
     228             : {
     229         515 :   PetscFunctionBegin;
     230         515 :   PetscCall(NEPInitializePackage());
     231         515 :   PetscCall(PetscFunctionListAdd(&NEPList,name,function));
     232         515 :   PetscFunctionReturn(PETSC_SUCCESS);
     233             : }
     234             : 
     235             : /*@C
     236             :    NEPMonitorRegister - Adds NEP monitor routine.
     237             : 
     238             :    Not Collective
     239             : 
     240             :    Input Parameters:
     241             : +  name    - name of a new monitor routine
     242             : .  vtype   - a PetscViewerType for the output
     243             : .  format  - a PetscViewerFormat for the output
     244             : .  monitor - monitor routine
     245             : .  create  - creation routine, or NULL
     246             : -  destroy - destruction routine, or NULL
     247             : 
     248             :    Notes:
     249             :    NEPMonitorRegister() may be called multiple times to add several user-defined monitors.
     250             : 
     251             :    Example Usage:
     252             : .vb
     253             :    NEPMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
     254             : .ve
     255             : 
     256             :    Then, your monitor can be chosen with the procedural interface via
     257             : $      NEPMonitorSetFromOptions(nep,"-nep_monitor_my_monitor","my_monitor",NULL)
     258             :    or at runtime via the option
     259             : $      -nep_monitor_my_monitor
     260             : 
     261             :    Level: advanced
     262             : 
     263             : .seealso: NEPMonitorRegisterAll()
     264             : @*/
     265         618 : PetscErrorCode NEPMonitorRegister(const char name[],PetscViewerType vtype,PetscViewerFormat format,PetscErrorCode (*monitor)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode (*create)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode (*destroy)(PetscViewerAndFormat**))
     266             : {
     267         618 :   char           key[PETSC_MAX_PATH_LEN];
     268             : 
     269         618 :   PetscFunctionBegin;
     270         618 :   PetscCall(NEPInitializePackage());
     271         618 :   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
     272         618 :   PetscCall(PetscFunctionListAdd(&NEPMonitorList,key,monitor));
     273         618 :   if (create)  PetscCall(PetscFunctionListAdd(&NEPMonitorCreateList,key,create));
     274         618 :   if (destroy) PetscCall(PetscFunctionListAdd(&NEPMonitorDestroyList,key,destroy));
     275         618 :   PetscFunctionReturn(PETSC_SUCCESS);
     276             : }
     277             : 
     278             : /*
     279             :    NEPReset_Problem - Destroys the problem matrices.
     280             : */
     281         192 : PetscErrorCode NEPReset_Problem(NEP nep)
     282             : {
     283         192 :   PetscInt       i;
     284             : 
     285         192 :   PetscFunctionBegin;
     286         192 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     287         192 :   PetscCall(MatDestroy(&nep->function));
     288         192 :   PetscCall(MatDestroy(&nep->function_pre));
     289         192 :   PetscCall(MatDestroy(&nep->jacobian));
     290         192 :   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
     291          85 :     PetscCall(MatDestroyMatrices(nep->nt,&nep->A));
     292         332 :     for (i=0;i<nep->nt;i++) PetscCall(FNDestroy(&nep->f[i]));
     293          85 :     PetscCall(PetscFree(nep->f));
     294          85 :     PetscCall(PetscFree(nep->nrma));
     295          85 :     if (nep->P) PetscCall(MatDestroyMatrices(nep->nt,&nep->P));
     296          85 :     nep->nt = 0;
     297             :   }
     298         192 :   PetscFunctionReturn(PETSC_SUCCESS);
     299             : }
     300             : /*@
     301             :    NEPReset - Resets the NEP context to the initial state (prior to setup)
     302             :    and destroys any allocated Vecs and Mats.
     303             : 
     304             :    Collective
     305             : 
     306             :    Input Parameter:
     307             : .  nep - eigensolver context obtained from NEPCreate()
     308             : 
     309             :    Level: advanced
     310             : 
     311             : .seealso: NEPDestroy()
     312             : @*/
     313         121 : PetscErrorCode NEPReset(NEP nep)
     314             : {
     315         121 :   PetscFunctionBegin;
     316         121 :   if (nep) PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     317         121 :   if (!nep) PetscFunctionReturn(PETSC_SUCCESS);
     318         121 :   PetscTryTypeMethod(nep,reset);
     319         121 :   if (nep->refineksp) PetscCall(KSPReset(nep->refineksp));
     320         121 :   PetscCall(NEPReset_Problem(nep));
     321         121 :   PetscCall(BVDestroy(&nep->V));
     322         121 :   PetscCall(BVDestroy(&nep->W));
     323         121 :   PetscCall(VecDestroyVecs(nep->nwork,&nep->work));
     324         121 :   PetscCall(MatDestroy(&nep->resolvent));
     325         121 :   nep->nwork = 0;
     326         121 :   nep->state = NEP_STATE_INITIAL;
     327         121 :   PetscFunctionReturn(PETSC_SUCCESS);
     328             : }
     329             : 
     330             : /*@
     331             :    NEPDestroy - Destroys the NEP context.
     332             : 
     333             :    Collective
     334             : 
     335             :    Input Parameter:
     336             : .  nep - eigensolver context obtained from NEPCreate()
     337             : 
     338             :    Level: beginner
     339             : 
     340             : .seealso: NEPCreate(), NEPSetUp(), NEPSolve()
     341             : @*/
     342         102 : PetscErrorCode NEPDestroy(NEP *nep)
     343             : {
     344         102 :   PetscFunctionBegin;
     345         102 :   if (!*nep) PetscFunctionReturn(PETSC_SUCCESS);
     346         102 :   PetscValidHeaderSpecific(*nep,NEP_CLASSID,1);
     347         102 :   if (--((PetscObject)*nep)->refct > 0) { *nep = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
     348         102 :   PetscCall(NEPReset(*nep));
     349         102 :   PetscTryTypeMethod(*nep,destroy);
     350         102 :   if ((*nep)->eigr) PetscCall(PetscFree4((*nep)->eigr,(*nep)->eigi,(*nep)->errest,(*nep)->perm));
     351         102 :   PetscCall(RGDestroy(&(*nep)->rg));
     352         102 :   PetscCall(DSDestroy(&(*nep)->ds));
     353         102 :   PetscCall(KSPDestroy(&(*nep)->refineksp));
     354         102 :   PetscCall(PetscSubcommDestroy(&(*nep)->refinesubc));
     355         102 :   PetscCall(PetscFree((*nep)->sc));
     356             :   /* just in case the initial vectors have not been used */
     357         102 :   PetscCall(SlepcBasisDestroy_Private(&(*nep)->nini,&(*nep)->IS));
     358         102 :   if ((*nep)->convergeddestroy) PetscCall((*(*nep)->convergeddestroy)(&(*nep)->convergedctx));
     359         102 :   if ((*nep)->stoppingdestroy) PetscCall((*(*nep)->stoppingdestroy)(&(*nep)->stoppingctx));
     360         102 :   PetscCall(NEPMonitorCancel(*nep));
     361         102 :   PetscCall(PetscHeaderDestroy(nep));
     362         102 :   PetscFunctionReturn(PETSC_SUCCESS);
     363             : }
     364             : 
     365             : /*@
     366             :    NEPSetBV - Associates a basis vectors object to the nonlinear eigensolver.
     367             : 
     368             :    Collective
     369             : 
     370             :    Input Parameters:
     371             : +  nep - eigensolver context obtained from NEPCreate()
     372             : -  bv  - the basis vectors object
     373             : 
     374             :    Note:
     375             :    Use NEPGetBV() to retrieve the basis vectors context (for example,
     376             :    to free it at the end of the computations).
     377             : 
     378             :    Level: advanced
     379             : 
     380             : .seealso: NEPGetBV()
     381             : @*/
     382           0 : PetscErrorCode NEPSetBV(NEP nep,BV bv)
     383             : {
     384           0 :   PetscFunctionBegin;
     385           0 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     386           0 :   PetscValidHeaderSpecific(bv,BV_CLASSID,2);
     387           0 :   PetscCheckSameComm(nep,1,bv,2);
     388           0 :   PetscCall(PetscObjectReference((PetscObject)bv));
     389           0 :   PetscCall(BVDestroy(&nep->V));
     390           0 :   nep->V = bv;
     391           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     392             : }
     393             : 
     394             : /*@
     395             :    NEPGetBV - Obtain the basis vectors object associated to the nonlinear
     396             :    eigensolver object.
     397             : 
     398             :    Not Collective
     399             : 
     400             :    Input Parameters:
     401             : .  nep - eigensolver context obtained from NEPCreate()
     402             : 
     403             :    Output Parameter:
     404             : .  bv - basis vectors context
     405             : 
     406             :    Level: advanced
     407             : 
     408             : .seealso: NEPSetBV()
     409             : @*/
     410         122 : PetscErrorCode NEPGetBV(NEP nep,BV *bv)
     411             : {
     412         122 :   PetscFunctionBegin;
     413         122 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     414         122 :   PetscAssertPointer(bv,2);
     415         122 :   if (!nep->V) {
     416         121 :     PetscCall(BVCreate(PetscObjectComm((PetscObject)nep),&nep->V));
     417         121 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)nep->V,(PetscObject)nep,0));
     418         121 :     PetscCall(PetscObjectSetOptions((PetscObject)nep->V,((PetscObject)nep)->options));
     419             :   }
     420         122 :   *bv = nep->V;
     421         122 :   PetscFunctionReturn(PETSC_SUCCESS);
     422             : }
     423             : 
     424             : /*@
     425             :    NEPSetRG - Associates a region object to the nonlinear eigensolver.
     426             : 
     427             :    Collective
     428             : 
     429             :    Input Parameters:
     430             : +  nep - eigensolver context obtained from NEPCreate()
     431             : -  rg  - the region object
     432             : 
     433             :    Note:
     434             :    Use NEPGetRG() to retrieve the region context (for example,
     435             :    to free it at the end of the computations).
     436             : 
     437             :    Level: advanced
     438             : 
     439             : .seealso: NEPGetRG()
     440             : @*/
     441           1 : PetscErrorCode NEPSetRG(NEP nep,RG rg)
     442             : {
     443           1 :   PetscFunctionBegin;
     444           1 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     445           1 :   if (rg) {
     446           1 :     PetscValidHeaderSpecific(rg,RG_CLASSID,2);
     447           1 :     PetscCheckSameComm(nep,1,rg,2);
     448             :   }
     449           1 :   PetscCall(PetscObjectReference((PetscObject)rg));
     450           1 :   PetscCall(RGDestroy(&nep->rg));
     451           1 :   nep->rg = rg;
     452           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     453             : }
     454             : 
     455             : /*@
     456             :    NEPGetRG - Obtain the region object associated to the
     457             :    nonlinear eigensolver object.
     458             : 
     459             :    Not Collective
     460             : 
     461             :    Input Parameters:
     462             : .  nep - eigensolver context obtained from NEPCreate()
     463             : 
     464             :    Output Parameter:
     465             : .  rg - region context
     466             : 
     467             :    Level: advanced
     468             : 
     469             : .seealso: NEPSetRG()
     470             : @*/
     471         102 : PetscErrorCode NEPGetRG(NEP nep,RG *rg)
     472             : {
     473         102 :   PetscFunctionBegin;
     474         102 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     475         102 :   PetscAssertPointer(rg,2);
     476         102 :   if (!nep->rg) {
     477         101 :     PetscCall(RGCreate(PetscObjectComm((PetscObject)nep),&nep->rg));
     478         101 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)nep->rg,(PetscObject)nep,0));
     479         101 :     PetscCall(PetscObjectSetOptions((PetscObject)nep->rg,((PetscObject)nep)->options));
     480             :   }
     481         102 :   *rg = nep->rg;
     482         102 :   PetscFunctionReturn(PETSC_SUCCESS);
     483             : }
     484             : 
     485             : /*@
     486             :    NEPSetDS - Associates a direct solver object to the nonlinear eigensolver.
     487             : 
     488             :    Collective
     489             : 
     490             :    Input Parameters:
     491             : +  nep - eigensolver context obtained from NEPCreate()
     492             : -  ds  - the direct solver object
     493             : 
     494             :    Note:
     495             :    Use NEPGetDS() to retrieve the direct solver context (for example,
     496             :    to free it at the end of the computations).
     497             : 
     498             :    Level: advanced
     499             : 
     500             : .seealso: NEPGetDS()
     501             : @*/
     502           0 : PetscErrorCode NEPSetDS(NEP nep,DS ds)
     503             : {
     504           0 :   PetscFunctionBegin;
     505           0 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     506           0 :   PetscValidHeaderSpecific(ds,DS_CLASSID,2);
     507           0 :   PetscCheckSameComm(nep,1,ds,2);
     508           0 :   PetscCall(PetscObjectReference((PetscObject)ds));
     509           0 :   PetscCall(DSDestroy(&nep->ds));
     510           0 :   nep->ds = ds;
     511           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     512             : }
     513             : 
     514             : /*@
     515             :    NEPGetDS - Obtain the direct solver object associated to the
     516             :    nonlinear eigensolver object.
     517             : 
     518             :    Not Collective
     519             : 
     520             :    Input Parameters:
     521             : .  nep - eigensolver context obtained from NEPCreate()
     522             : 
     523             :    Output Parameter:
     524             : .  ds - direct solver context
     525             : 
     526             :    Level: advanced
     527             : 
     528             : .seealso: NEPSetDS()
     529             : @*/
     530          80 : PetscErrorCode NEPGetDS(NEP nep,DS *ds)
     531             : {
     532          80 :   PetscFunctionBegin;
     533          80 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     534          80 :   PetscAssertPointer(ds,2);
     535          80 :   if (!nep->ds) {
     536          79 :     PetscCall(DSCreate(PetscObjectComm((PetscObject)nep),&nep->ds));
     537          79 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)nep->ds,(PetscObject)nep,0));
     538          79 :     PetscCall(PetscObjectSetOptions((PetscObject)nep->ds,((PetscObject)nep)->options));
     539             :   }
     540          80 :   *ds = nep->ds;
     541          80 :   PetscFunctionReturn(PETSC_SUCCESS);
     542             : }
     543             : 
     544             : /*@
     545             :    NEPRefineGetKSP - Obtain the ksp object used by the eigensolver
     546             :    object in the refinement phase.
     547             : 
     548             :    Collective
     549             : 
     550             :    Input Parameters:
     551             : .  nep - eigensolver context obtained from NEPCreate()
     552             : 
     553             :    Output Parameter:
     554             : .  ksp - ksp context
     555             : 
     556             :    Level: advanced
     557             : 
     558             : .seealso: NEPSetRefine()
     559             : @*/
     560         106 : PetscErrorCode NEPRefineGetKSP(NEP nep,KSP *ksp)
     561             : {
     562         106 :   MPI_Comm       comm;
     563             : 
     564         106 :   PetscFunctionBegin;
     565         106 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     566         106 :   PetscAssertPointer(ksp,2);
     567         106 :   if (!nep->refineksp) {
     568         102 :     if (nep->npart>1) {
     569             :       /* Split in subcomunicators */
     570           8 :       PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&nep->refinesubc));
     571           8 :       PetscCall(PetscSubcommSetNumber(nep->refinesubc,nep->npart));
     572           8 :       PetscCall(PetscSubcommSetType(nep->refinesubc,PETSC_SUBCOMM_CONTIGUOUS));
     573           8 :       PetscCall(PetscSubcommGetChild(nep->refinesubc,&comm));
     574          94 :     } else PetscCall(PetscObjectGetComm((PetscObject)nep,&comm));
     575         102 :     PetscCall(KSPCreate(comm,&nep->refineksp));
     576         102 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)nep->refineksp,(PetscObject)nep,0));
     577         102 :     PetscCall(PetscObjectSetOptions((PetscObject)nep->refineksp,((PetscObject)nep)->options));
     578         102 :     PetscCall(KSPSetOptionsPrefix(*ksp,((PetscObject)nep)->prefix));
     579         102 :     PetscCall(KSPAppendOptionsPrefix(*ksp,"nep_refine_"));
     580         203 :     PetscCall(KSPSetTolerances(nep->refineksp,SlepcDefaultTol(nep->rtol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
     581             :   }
     582         106 :   *ksp = nep->refineksp;
     583         106 :   PetscFunctionReturn(PETSC_SUCCESS);
     584             : }
     585             : 
     586             : /*@
     587             :    NEPSetTarget - Sets the value of the target.
     588             : 
     589             :    Logically Collective
     590             : 
     591             :    Input Parameters:
     592             : +  nep    - eigensolver context
     593             : -  target - the value of the target
     594             : 
     595             :    Options Database Key:
     596             : .  -nep_target <scalar> - the value of the target
     597             : 
     598             :    Notes:
     599             :    The target is a scalar value used to determine the portion of the spectrum
     600             :    of interest. It is used in combination with NEPSetWhichEigenpairs().
     601             : 
     602             :    In the case of complex scalars, a complex value can be provided in the
     603             :    command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
     604             :    -nep_target 1.0+2.0i
     605             : 
     606             :    Level: intermediate
     607             : 
     608             : .seealso: NEPGetTarget(), NEPSetWhichEigenpairs()
     609             : @*/
     610          89 : PetscErrorCode NEPSetTarget(NEP nep,PetscScalar target)
     611             : {
     612          89 :   PetscFunctionBegin;
     613          89 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     614         267 :   PetscValidLogicalCollectiveScalar(nep,target,2);
     615          89 :   nep->target = target;
     616          89 :   PetscFunctionReturn(PETSC_SUCCESS);
     617             : }
     618             : 
     619             : /*@
     620             :    NEPGetTarget - Gets the value of the target.
     621             : 
     622             :    Not Collective
     623             : 
     624             :    Input Parameter:
     625             : .  nep - eigensolver context
     626             : 
     627             :    Output Parameter:
     628             : .  target - the value of the target
     629             : 
     630             :    Note:
     631             :    If the target was not set by the user, then zero is returned.
     632             : 
     633             :    Level: intermediate
     634             : 
     635             : .seealso: NEPSetTarget()
     636             : @*/
     637           2 : PetscErrorCode NEPGetTarget(NEP nep,PetscScalar* target)
     638             : {
     639           2 :   PetscFunctionBegin;
     640           2 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     641           2 :   PetscAssertPointer(target,2);
     642           2 :   *target = nep->target;
     643           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     644             : }
     645             : 
     646             : /*@C
     647             :    NEPSetFunction - Sets the function to compute the nonlinear Function T(lambda)
     648             :    as well as the location to store the matrix.
     649             : 
     650             :    Collective
     651             : 
     652             :    Input Parameters:
     653             : +  nep - the NEP context
     654             : .  A   - Function matrix
     655             : .  B   - preconditioner matrix (usually same as A)
     656             : .  fun - Function evaluation routine (if NULL then NEP retains any
     657             :          previously set value), see NEPFunctionFn for the calling sequence
     658             : -  ctx - [optional] user-defined context for private data for the Function
     659             :          evaluation routine (may be NULL) (if NULL then NEP retains any
     660             :          previously set value)
     661             : 
     662             :    Level: beginner
     663             : 
     664             : .seealso: NEPGetFunction(), NEPSetJacobian()
     665             : @*/
     666          36 : PetscErrorCode NEPSetFunction(NEP nep,Mat A,Mat B,NEPFunctionFn *fun,void *ctx)
     667             : {
     668          36 :   PetscFunctionBegin;
     669          36 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     670          36 :   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
     671          36 :   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
     672          36 :   if (A) PetscCheckSameComm(nep,1,A,2);
     673          36 :   if (B) PetscCheckSameComm(nep,1,B,3);
     674             : 
     675          36 :   if (nep->state) PetscCall(NEPReset(nep));
     676          31 :   else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPReset_Problem(nep));
     677             : 
     678          36 :   if (fun) nep->computefunction = fun;
     679          36 :   if (ctx) nep->functionctx     = ctx;
     680          36 :   if (A) {
     681          36 :     PetscCall(PetscObjectReference((PetscObject)A));
     682          36 :     PetscCall(MatDestroy(&nep->function));
     683          36 :     nep->function = A;
     684             :   }
     685          36 :   if (B) {
     686          36 :     PetscCall(PetscObjectReference((PetscObject)B));
     687          36 :     PetscCall(MatDestroy(&nep->function_pre));
     688          36 :     nep->function_pre = B;
     689             :   }
     690          36 :   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
     691          36 :   nep->state = NEP_STATE_INITIAL;
     692          36 :   PetscFunctionReturn(PETSC_SUCCESS);
     693             : }
     694             : 
     695             : /*@C
     696             :    NEPGetFunction - Returns the Function matrix and optionally the user
     697             :    provided context for evaluating the Function.
     698             : 
     699             :    Not Collective
     700             : 
     701             :    Input Parameter:
     702             : .  nep - the nonlinear eigensolver context
     703             : 
     704             :    Output Parameters:
     705             : +  A   - location to stash Function matrix (or NULL)
     706             : .  B   - location to stash preconditioner matrix (or NULL)
     707             : .  fun - location to put Function function (or NULL)
     708             : -  ctx - location to stash Function context (or NULL)
     709             : 
     710             :    Level: advanced
     711             : 
     712             : .seealso: NEPSetFunction()
     713             : @*/
     714          72 : PetscErrorCode NEPGetFunction(NEP nep,Mat *A,Mat *B,NEPFunctionFn **fun,void **ctx)
     715             : {
     716          72 :   PetscFunctionBegin;
     717          72 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     718          72 :   NEPCheckCallback(nep,1);
     719          72 :   if (A)   *A   = nep->function;
     720          72 :   if (B)   *B   = nep->function_pre;
     721          72 :   if (fun) *fun = nep->computefunction;
     722          72 :   if (ctx) *ctx = nep->functionctx;
     723          72 :   PetscFunctionReturn(PETSC_SUCCESS);
     724             : }
     725             : 
     726             : /*@C
     727             :    NEPSetJacobian - Sets the function to compute the Jacobian T'(lambda) as well
     728             :    as the location to store the matrix.
     729             : 
     730             :    Collective
     731             : 
     732             :    Input Parameters:
     733             : +  nep - the NEP context
     734             : .  A   - Jacobian matrix
     735             : .  jac - Jacobian evaluation routine (if NULL then NEP retains any
     736             :          previously set value), see NEPJacobianFn for the calling sequence
     737             : -  ctx - [optional] user-defined context for private data for the Jacobian
     738             :          evaluation routine (may be NULL) (if NULL then NEP retains any
     739             :          previously set value)
     740             : 
     741             :    Level: beginner
     742             : 
     743             : .seealso: NEPSetFunction(), NEPGetJacobian()
     744             : @*/
     745          33 : PetscErrorCode NEPSetJacobian(NEP nep,Mat A,NEPJacobianFn *jac,void *ctx)
     746             : {
     747          33 :   PetscFunctionBegin;
     748          33 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     749          33 :   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
     750          33 :   if (A) PetscCheckSameComm(nep,1,A,2);
     751             : 
     752          33 :   if (nep->state) PetscCall(NEPReset(nep));
     753          33 :   else if (nep->fui && nep->fui!=NEP_USER_INTERFACE_CALLBACK) PetscCall(NEPReset_Problem(nep));
     754             : 
     755          33 :   if (jac) nep->computejacobian = jac;
     756          33 :   if (ctx) nep->jacobianctx     = ctx;
     757          33 :   if (A) {
     758          33 :     PetscCall(PetscObjectReference((PetscObject)A));
     759          33 :     PetscCall(MatDestroy(&nep->jacobian));
     760          33 :     nep->jacobian = A;
     761             :   }
     762          33 :   nep->fui   = NEP_USER_INTERFACE_CALLBACK;
     763          33 :   nep->state = NEP_STATE_INITIAL;
     764          33 :   PetscFunctionReturn(PETSC_SUCCESS);
     765             : }
     766             : 
     767             : /*@C
     768             :    NEPGetJacobian - Returns the Jacobian matrix and optionally the user
     769             :    provided routine and context for evaluating the Jacobian.
     770             : 
     771             :    Not Collective
     772             : 
     773             :    Input Parameter:
     774             : .  nep - the nonlinear eigensolver context
     775             : 
     776             :    Output Parameters:
     777             : +  A   - location to stash Jacobian matrix (or NULL)
     778             : .  jac - location to put Jacobian function (or NULL)
     779             : -  ctx - location to stash Jacobian context (or NULL)
     780             : 
     781             :    Level: advanced
     782             : 
     783             : .seealso: NEPSetJacobian()
     784             : @*/
     785           0 : PetscErrorCode NEPGetJacobian(NEP nep,Mat *A,NEPJacobianFn **jac,void **ctx)
     786             : {
     787           0 :   PetscFunctionBegin;
     788           0 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     789           0 :   NEPCheckCallback(nep,1);
     790           0 :   if (A)   *A   = nep->jacobian;
     791           0 :   if (jac) *jac = nep->computejacobian;
     792           0 :   if (ctx) *ctx = nep->jacobianctx;
     793           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     794             : }
     795             : 
     796             : /*@
     797             :    NEPSetSplitOperator - Sets the operator of the nonlinear eigenvalue problem
     798             :    in split form.
     799             : 
     800             :    Collective
     801             : 
     802             :    Input Parameters:
     803             : +  nep - the nonlinear eigensolver context
     804             : .  nt  - number of terms in the split form
     805             : .  A   - array of matrices
     806             : .  f   - array of functions
     807             : -  str - structure flag for matrices
     808             : 
     809             :    Notes:
     810             :    The nonlinear operator is written as T(lambda) = sum_i A_i*f_i(lambda),
     811             :    for i=1,...,n. The derivative T'(lambda) can be obtained using the
     812             :    derivatives of f_i.
     813             : 
     814             :    The structure flag provides information about A_i's nonzero pattern
     815             :    (see MatStructure enum). If all matrices have the same pattern, then
     816             :    use SAME_NONZERO_PATTERN. If the patterns are different but contained
     817             :    in the pattern of the first one, then use SUBSET_NONZERO_PATTERN. If
     818             :    patterns are known to be different, use DIFFERENT_NONZERO_PATTERN.
     819             :    If set to UNKNOWN_NONZERO_PATTERN, the patterns will be compared to
     820             :    determine if they are equal.
     821             : 
     822             :    This function must be called before NEPSetUp(). If it is called again
     823             :    after NEPSetUp() then the NEP object is reset.
     824             : 
     825             :    Level: beginner
     826             : 
     827             : .seealso: NEPGetSplitOperatorTerm(), NEPGetSplitOperatorInfo(), NEPSetSplitPreconditioner()
     828             : @*/
     829          85 : PetscErrorCode NEPSetSplitOperator(NEP nep,PetscInt nt,Mat A[],FN f[],MatStructure str)
     830             : {
     831          85 :   PetscInt       i,n=0,m,m0=0,mloc,nloc,mloc0=0;
     832             : 
     833          85 :   PetscFunctionBegin;
     834          85 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     835         255 :   PetscValidLogicalCollectiveInt(nep,nt,2);
     836          85 :   PetscCheck(nt>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more terms, you have %" PetscInt_FMT,nt);
     837          85 :   PetscAssertPointer(A,3);
     838          85 :   PetscAssertPointer(f,4);
     839         255 :   PetscValidLogicalCollectiveEnum(nep,str,5);
     840             : 
     841         332 :   for (i=0;i<nt;i++) {
     842         247 :     PetscValidHeaderSpecific(A[i],MAT_CLASSID,3);
     843         247 :     PetscCheckSameComm(nep,1,A[i],3);
     844         247 :     PetscValidHeaderSpecific(f[i],FN_CLASSID,4);
     845         247 :     PetscCheckSameComm(nep,1,f[i],4);
     846         247 :     PetscCall(MatGetSize(A[i],&m,&n));
     847         247 :     PetscCall(MatGetLocalSize(A[i],&mloc,&nloc));
     848         247 :     PetscCheck(m==n,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"A[%" PetscInt_FMT "] is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",i,m,n);
     849         247 :     PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"A[%" PetscInt_FMT "] does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,nloc);
     850         247 :     if (!i) { m0 = m; mloc0 = mloc; }
     851         247 :     PetscCheck(m==m0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Dimensions of A[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,m,m0);
     852         247 :     PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Local dimensions of A[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,mloc0);
     853         247 :     PetscCall(PetscObjectReference((PetscObject)A[i]));
     854         247 :     PetscCall(PetscObjectReference((PetscObject)f[i]));
     855             :   }
     856             : 
     857          85 :   if (nep->state && (n!=nep->n || nloc!=nep->nloc)) PetscCall(NEPReset(nep));
     858          71 :   else PetscCall(NEPReset_Problem(nep));
     859             : 
     860             :   /* allocate space and copy matrices and functions */
     861          85 :   PetscCall(PetscMalloc1(nt,&nep->A));
     862         332 :   for (i=0;i<nt;i++) nep->A[i] = A[i];
     863          85 :   PetscCall(PetscMalloc1(nt,&nep->f));
     864         332 :   for (i=0;i<nt;i++) nep->f[i] = f[i];
     865          85 :   PetscCall(PetscCalloc1(nt,&nep->nrma));
     866          85 :   nep->nt    = nt;
     867          85 :   nep->mstr  = str;
     868          85 :   nep->fui   = NEP_USER_INTERFACE_SPLIT;
     869          85 :   nep->state = NEP_STATE_INITIAL;
     870          85 :   PetscFunctionReturn(PETSC_SUCCESS);
     871             : }
     872             : 
     873             : /*@
     874             :    NEPGetSplitOperatorTerm - Gets the matrices and functions associated with
     875             :    the nonlinear operator in split form.
     876             : 
     877             :    Not Collective
     878             : 
     879             :    Input Parameters:
     880             : +  nep - the nonlinear eigensolver context
     881             : -  k   - the index of the requested term (starting in 0)
     882             : 
     883             :    Output Parameters:
     884             : +  A - the matrix of the requested term
     885             : -  f - the function of the requested term
     886             : 
     887             :    Level: intermediate
     888             : 
     889             : .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorInfo()
     890             : @*/
     891          25 : PetscErrorCode NEPGetSplitOperatorTerm(NEP nep,PetscInt k,Mat *A,FN *f)
     892             : {
     893          25 :   PetscFunctionBegin;
     894          25 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     895          75 :   PetscValidLogicalCollectiveInt(nep,k,2);
     896          25 :   NEPCheckSplit(nep,1);
     897          25 :   PetscCheck(k>=0 && k<nep->nt,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,nep->nt-1);
     898          25 :   if (A) *A = nep->A[k];
     899          25 :   if (f) *f = nep->f[k];
     900          25 :   PetscFunctionReturn(PETSC_SUCCESS);
     901             : }
     902             : 
     903             : /*@
     904             :    NEPGetSplitOperatorInfo - Returns the number of terms of the split form of
     905             :    the nonlinear operator, as well as the structure flag for matrices.
     906             : 
     907             :    Not Collective
     908             : 
     909             :    Input Parameter:
     910             : .  nep - the nonlinear eigensolver context
     911             : 
     912             :    Output Parameters:
     913             : +  n   - the number of terms passed in NEPSetSplitOperator()
     914             : -  str - the matrix structure flag passed in NEPSetSplitOperator()
     915             : 
     916             :    Level: intermediate
     917             : 
     918             : .seealso: NEPSetSplitOperator(), NEPGetSplitOperatorTerm()
     919             : @*/
     920           9 : PetscErrorCode NEPGetSplitOperatorInfo(NEP nep,PetscInt *n,MatStructure *str)
     921             : {
     922           9 :   PetscFunctionBegin;
     923           9 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     924           9 :   NEPCheckSplit(nep,1);
     925           9 :   if (n)   *n = nep->nt;
     926           9 :   if (str) *str = nep->mstr;
     927           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     928             : }
     929             : 
     930             : /*@
     931             :    NEPSetSplitPreconditioner - Sets an operator in split form from which
     932             :    to build the preconditioner to be used when solving the nonlinear
     933             :    eigenvalue problem in split form.
     934             : 
     935             :    Collective
     936             : 
     937             :    Input Parameters:
     938             : +  nep  - the nonlinear eigensolver context
     939             : .  ntp  - number of terms in the split preconditioner
     940             : .  P    - array of matrices
     941             : -  strp - structure flag for matrices
     942             : 
     943             :    Notes:
     944             :    The matrix for the preconditioner is expressed as P(lambda) =
     945             :    sum_i P_i*f_i(lambda), for i=1,...,n, where the f_i functions
     946             :    are the same as in NEPSetSplitOperator(). It is not necessary to call
     947             :    this function. If it is not invoked, then the preconditioner is
     948             :    built from T(lambda), i.e., both matrices and functions passed in
     949             :    NEPSetSplitOperator().
     950             : 
     951             :    The structure flag provides information about P_i's nonzero pattern
     952             :    in the same way as in NEPSetSplitOperator().
     953             : 
     954             :    If the functions defining the preconditioner operator were different
     955             :    from the ones given in NEPSetSplitOperator(), then the split form
     956             :    cannot be used. Use the callback interface instead.
     957             : 
     958             :    Use ntp=0 to reset a previously set split preconditioner.
     959             : 
     960             :    Level: advanced
     961             : 
     962             : .seealso: NEPGetSplitPreconditionerTerm(), NEPGetSplitPreconditionerInfo(), NEPSetSplitOperator()
     963             : @*/
     964           4 : PetscErrorCode NEPSetSplitPreconditioner(NEP nep,PetscInt ntp,Mat P[],MatStructure strp)
     965             : {
     966           4 :   PetscInt       i,n=0,m,m0=0,mloc,nloc,mloc0=0;
     967             : 
     968           4 :   PetscFunctionBegin;
     969           4 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
     970          12 :   PetscValidLogicalCollectiveInt(nep,ntp,2);
     971           4 :   PetscCheck(ntp>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Negative value of ntp = %" PetscInt_FMT,ntp);
     972           4 :   PetscCheck(nep->fui==NEP_USER_INTERFACE_SPLIT,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetSplitOperator first");
     973           4 :   PetscCheck(ntp==0 || nep->nt==ntp,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The number of terms must be the same as in NEPSetSplitOperator()");
     974           4 :   if (ntp) PetscAssertPointer(P,3);
     975          12 :   PetscValidLogicalCollectiveEnum(nep,strp,4);
     976             : 
     977          16 :   for (i=0;i<ntp;i++) {
     978          12 :     PetscValidHeaderSpecific(P[i],MAT_CLASSID,3);
     979          12 :     PetscCheckSameComm(nep,1,P[i],3);
     980          12 :     PetscCall(MatGetSize(P[i],&m,&n));
     981          12 :     PetscCall(MatGetLocalSize(P[i],&mloc,&nloc));
     982          12 :     PetscCheck(m==n,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"P[%" PetscInt_FMT "] is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",i,m,n);
     983          12 :     PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"P[%" PetscInt_FMT "] does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,nloc);
     984          12 :     if (!i) { m0 = m; mloc0 = mloc; }
     985          12 :     PetscCheck(m==m0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Dimensions of P[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,m,m0);
     986          12 :     PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_INCOMP,"Local dimensions of P[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,mloc0);
     987          12 :     PetscCall(PetscObjectReference((PetscObject)P[i]));
     988             :   }
     989             : 
     990           4 :   PetscCheck(!nep->state,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"To call this function after NEPSetUp(), you must call NEPSetSplitOperator() again");
     991           4 :   if (nep->P) PetscCall(MatDestroyMatrices(nep->nt,&nep->P));
     992             : 
     993             :   /* allocate space and copy matrices */
     994           4 :   if (ntp) {
     995           4 :     PetscCall(PetscMalloc1(ntp,&nep->P));
     996          16 :     for (i=0;i<ntp;i++) nep->P[i] = P[i];
     997             :   }
     998           4 :   nep->mstrp = strp;
     999           4 :   nep->state = NEP_STATE_INITIAL;
    1000           4 :   PetscFunctionReturn(PETSC_SUCCESS);
    1001             : }
    1002             : 
    1003             : /*@
    1004             :    NEPGetSplitPreconditionerTerm - Gets the matrices associated with
    1005             :    the split preconditioner.
    1006             : 
    1007             :    Not Collective
    1008             : 
    1009             :    Input Parameters:
    1010             : +  nep - the nonlinear eigensolver context
    1011             : -  k   - the index of the requested term (starting in 0)
    1012             : 
    1013             :    Output Parameter:
    1014             : .  P  - the matrix of the requested term
    1015             : 
    1016             :    Level: advanced
    1017             : 
    1018             : .seealso: NEPSetSplitPreconditioner(), NEPGetSplitPreconditionerInfo()
    1019             : @*/
    1020           4 : PetscErrorCode NEPGetSplitPreconditionerTerm(NEP nep,PetscInt k,Mat *P)
    1021             : {
    1022           4 :   PetscFunctionBegin;
    1023           4 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
    1024          12 :   PetscValidLogicalCollectiveInt(nep,k,2);
    1025           4 :   PetscAssertPointer(P,3);
    1026           4 :   NEPCheckSplit(nep,1);
    1027           4 :   PetscCheck(k>=0 && k<nep->nt,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,nep->nt-1);
    1028           4 :   PetscCheck(nep->P,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"You have not called NEPSetSplitPreconditioner()");
    1029           4 :   *P = nep->P[k];
    1030           4 :   PetscFunctionReturn(PETSC_SUCCESS);
    1031             : }
    1032             : 
    1033             : /*@
    1034             :    NEPGetSplitPreconditionerInfo - Returns the number of terms of the split
    1035             :    preconditioner, as well as the structure flag for matrices.
    1036             : 
    1037             :    Not Collective
    1038             : 
    1039             :    Input Parameter:
    1040             : .  nep - the nonlinear eigensolver context
    1041             : 
    1042             :    Output Parameters:
    1043             : +  n    - the number of terms passed in NEPSetSplitPreconditioner()
    1044             : -  strp - the matrix structure flag passed in NEPSetSplitPreconditioner()
    1045             : 
    1046             :    Level: advanced
    1047             : 
    1048             : .seealso: NEPSetSplitPreconditioner(), NEPGetSplitPreconditionerTerm()
    1049             : @*/
    1050           4 : PetscErrorCode NEPGetSplitPreconditionerInfo(NEP nep,PetscInt *n,MatStructure *strp)
    1051             : {
    1052           4 :   PetscFunctionBegin;
    1053           4 :   PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
    1054           4 :   NEPCheckSplit(nep,1);
    1055           4 :   if (n)    *n    = nep->P? nep->nt: 0;
    1056           4 :   if (strp) *strp = nep->mstrp;
    1057           4 :   PetscFunctionReturn(PETSC_SUCCESS);
    1058             : }

Generated by: LCOV version 1.14