1*f17f7ee2SSatish Balay #include <petsc/private/matimpl.h> 2*f17f7ee2SSatish Balay #include <petsc/private/vecimpl.h> 3*f17f7ee2SSatish Balay #include <petscsf.h> 4*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 5*f17f7ee2SSatish Balay #include <thrust/for_each.h> 6*f17f7ee2SSatish Balay #include <thrust/device_vector.h> 7*f17f7ee2SSatish Balay #include <thrust/execution_policy.h> 8*f17f7ee2SSatish Balay #endif 9*f17f7ee2SSatish Balay 10*f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode PetscSFGetVectorSF(PetscSF sf, PetscInt nv, PetscInt ldr, PetscInt ldl, PetscSF *vsf) 11*f17f7ee2SSatish Balay { 12*f17f7ee2SSatish Balay PetscSF rankssf; 13*f17f7ee2SSatish Balay const PetscSFNode *iremote; 14*f17f7ee2SSatish Balay PetscSFNode *viremote,*rremotes; 15*f17f7ee2SSatish Balay const PetscInt *ilocal; 16*f17f7ee2SSatish Balay PetscInt *vilocal = NULL,*ldrs; 17*f17f7ee2SSatish Balay const PetscMPIInt *ranks; 18*f17f7ee2SSatish Balay PetscMPIInt *sranks; 19*f17f7ee2SSatish Balay PetscInt nranks,nr,nl,vnr,vnl,i,v,j,maxl; 20*f17f7ee2SSatish Balay MPI_Comm comm; 21*f17f7ee2SSatish Balay 22*f17f7ee2SSatish Balay PetscFunctionBegin; 23*f17f7ee2SSatish Balay PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 24*f17f7ee2SSatish Balay PetscValidLogicalCollectiveInt(sf,nv,2); 25*f17f7ee2SSatish Balay PetscValidPointer(vsf,5); 26*f17f7ee2SSatish Balay if (nv == 1) { 27*f17f7ee2SSatish Balay PetscCall(PetscObjectReference((PetscObject)sf)); 28*f17f7ee2SSatish Balay *vsf = sf; 29*f17f7ee2SSatish Balay PetscFunctionReturn(0); 30*f17f7ee2SSatish Balay } 31*f17f7ee2SSatish Balay PetscCall(PetscObjectGetComm((PetscObject)sf,&comm)); 32*f17f7ee2SSatish Balay PetscCall(PetscSFGetGraph(sf,&nr,&nl,&ilocal,&iremote)); 33*f17f7ee2SSatish Balay PetscCall(PetscSFGetLeafRange(sf,NULL,&maxl)); 34*f17f7ee2SSatish Balay maxl += 1; 35*f17f7ee2SSatish Balay if (ldl == PETSC_DECIDE) ldl = maxl; 36*f17f7ee2SSatish Balay if (ldr == PETSC_DECIDE) ldr = nr; 37*f17f7ee2SSatish Balay PetscCheck(ldr >= nr,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid leading dimension %" PetscInt_FMT " < %" PetscInt_FMT,ldr,nr); 38*f17f7ee2SSatish Balay PetscCheck(ldl >= maxl,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid leading dimension %" PetscInt_FMT " < %" PetscInt_FMT,ldl,maxl); 39*f17f7ee2SSatish Balay vnr = nr*nv; 40*f17f7ee2SSatish Balay vnl = nl*nv; 41*f17f7ee2SSatish Balay PetscCall(PetscMalloc1(vnl,&viremote)); 42*f17f7ee2SSatish Balay if (ilocal) PetscCall(PetscMalloc1(vnl,&vilocal)); 43*f17f7ee2SSatish Balay 44*f17f7ee2SSatish Balay /* TODO: Should this special SF be available, e.g. 45*f17f7ee2SSatish Balay PetscSFGetRanksSF or similar? */ 46*f17f7ee2SSatish Balay PetscCall(PetscSFGetRootRanks(sf,&nranks,&ranks,NULL,NULL,NULL)); 47*f17f7ee2SSatish Balay PetscCall(PetscMalloc1(nranks,&sranks)); 48*f17f7ee2SSatish Balay PetscCall(PetscArraycpy(sranks,ranks,nranks)); 49*f17f7ee2SSatish Balay PetscCall(PetscSortMPIInt(nranks,sranks)); 50*f17f7ee2SSatish Balay PetscCall(PetscMalloc1(nranks,&rremotes)); 51*f17f7ee2SSatish Balay for (i=0;i<nranks;i++) { 52*f17f7ee2SSatish Balay rremotes[i].rank = sranks[i]; 53*f17f7ee2SSatish Balay rremotes[i].index = 0; 54*f17f7ee2SSatish Balay } 55*f17f7ee2SSatish Balay PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&rankssf)); 56*f17f7ee2SSatish Balay PetscCall(PetscSFSetGraph(rankssf,1,nranks,NULL,PETSC_OWN_POINTER,rremotes,PETSC_OWN_POINTER)); 57*f17f7ee2SSatish Balay PetscCall(PetscMalloc1(nranks,&ldrs)); 58*f17f7ee2SSatish Balay PetscCall(PetscSFBcastBegin(rankssf,MPIU_INT,&ldr,ldrs,MPI_REPLACE)); 59*f17f7ee2SSatish Balay PetscCall(PetscSFBcastEnd(rankssf,MPIU_INT,&ldr,ldrs,MPI_REPLACE)); 60*f17f7ee2SSatish Balay PetscCall(PetscSFDestroy(&rankssf)); 61*f17f7ee2SSatish Balay 62*f17f7ee2SSatish Balay j = -1; 63*f17f7ee2SSatish Balay for (i=0;i<nl;i++) { 64*f17f7ee2SSatish Balay const PetscInt r = iremote[i].rank; 65*f17f7ee2SSatish Balay const PetscInt ii = iremote[i].index; 66*f17f7ee2SSatish Balay 67*f17f7ee2SSatish Balay if (j < 0 || sranks[j] != r) { 68*f17f7ee2SSatish Balay PetscCall(PetscFindMPIInt(r,nranks,sranks,&j)); 69*f17f7ee2SSatish Balay } 70*f17f7ee2SSatish Balay PetscCheck(j >= 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unable to locate neighbor rank %" PetscInt_FMT,r); 71*f17f7ee2SSatish Balay for (v=0;v<nv;v++) { 72*f17f7ee2SSatish Balay viremote[v*nl + i].rank = r; 73*f17f7ee2SSatish Balay viremote[v*nl + i].index = v*ldrs[j] + ii; 74*f17f7ee2SSatish Balay if (ilocal) vilocal[v*nl + i] = v*ldl + ilocal[i]; 75*f17f7ee2SSatish Balay } 76*f17f7ee2SSatish Balay } 77*f17f7ee2SSatish Balay PetscCall(PetscFree(sranks)); 78*f17f7ee2SSatish Balay PetscCall(PetscFree(ldrs)); 79*f17f7ee2SSatish Balay PetscCall(PetscSFCreate(comm,vsf)); 80*f17f7ee2SSatish Balay PetscCall(PetscSFSetGraph(*vsf,vnr,vnl,vilocal,PETSC_OWN_POINTER,viremote,PETSC_OWN_POINTER)); 81*f17f7ee2SSatish Balay PetscFunctionReturn(0); 82*f17f7ee2SSatish Balay } 83*f17f7ee2SSatish Balay 84*f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode MatDenseGetH2OpusVectorSF(Mat A, PetscSF h2sf, PetscSF *osf) 85*f17f7ee2SSatish Balay { 86*f17f7ee2SSatish Balay PetscSF asf; 87*f17f7ee2SSatish Balay 88*f17f7ee2SSatish Balay PetscFunctionBegin; 89*f17f7ee2SSatish Balay PetscValidHeaderSpecific(A,MAT_CLASSID,1); 90*f17f7ee2SSatish Balay PetscValidHeaderSpecific(h2sf,PETSCSF_CLASSID,2); 91*f17f7ee2SSatish Balay PetscValidPointer(osf,3); 92*f17f7ee2SSatish Balay PetscCall(PetscObjectQuery((PetscObject)A,"_math2opus_vectorsf",(PetscObject*)&asf)); 93*f17f7ee2SSatish Balay if (!asf) { 94*f17f7ee2SSatish Balay PetscInt lda; 95*f17f7ee2SSatish Balay 96*f17f7ee2SSatish Balay PetscCall(MatDenseGetLDA(A,&lda)); 97*f17f7ee2SSatish Balay PetscCall(PetscSFGetVectorSF(h2sf,A->cmap->N,lda,PETSC_DECIDE,&asf)); 98*f17f7ee2SSatish Balay PetscCall(PetscObjectCompose((PetscObject)A,"_math2opus_vectorsf",(PetscObject)asf)); 99*f17f7ee2SSatish Balay PetscCall(PetscObjectDereference((PetscObject)asf)); 100*f17f7ee2SSatish Balay } 101*f17f7ee2SSatish Balay *osf = asf; 102*f17f7ee2SSatish Balay PetscFunctionReturn(0); 103*f17f7ee2SSatish Balay } 104*f17f7ee2SSatish Balay 105*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 106*f17f7ee2SSatish Balay struct SignVector_Functor 107*f17f7ee2SSatish Balay { 108*f17f7ee2SSatish Balay const PetscScalar *v; 109*f17f7ee2SSatish Balay PetscScalar *s; 110*f17f7ee2SSatish Balay SignVector_Functor(const PetscScalar *_v, PetscScalar *_s) : v(_v), s(_s) {} 111*f17f7ee2SSatish Balay 112*f17f7ee2SSatish Balay __host__ __device__ void operator()(PetscInt i) 113*f17f7ee2SSatish Balay { 114*f17f7ee2SSatish Balay s[i] = (v[i] < 0 ? -1 : 1); 115*f17f7ee2SSatish Balay } 116*f17f7ee2SSatish Balay }; 117*f17f7ee2SSatish Balay #endif 118*f17f7ee2SSatish Balay 119*f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSign(Vec v, Vec s) 120*f17f7ee2SSatish Balay { 121*f17f7ee2SSatish Balay const PetscScalar *av; 122*f17f7ee2SSatish Balay PetscScalar *as; 123*f17f7ee2SSatish Balay PetscInt i,n; 124*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 125*f17f7ee2SSatish Balay PetscBool viscuda,siscuda; 126*f17f7ee2SSatish Balay #endif 127*f17f7ee2SSatish Balay 128*f17f7ee2SSatish Balay PetscFunctionBegin; 129*f17f7ee2SSatish Balay PetscValidHeaderSpecific(v,VEC_CLASSID,1); 130*f17f7ee2SSatish Balay PetscValidHeaderSpecific(s,VEC_CLASSID,2); 131*f17f7ee2SSatish Balay PetscCall(VecGetLocalSize(s,&n)); 132*f17f7ee2SSatish Balay PetscCall(VecGetLocalSize(v,&i)); 133*f17f7ee2SSatish Balay PetscCheck(i == n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Invalid local sizes %" PetscInt_FMT " != %" PetscInt_FMT,i,n); 134*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 135*f17f7ee2SSatish Balay PetscCall(PetscObjectTypeCompareAny((PetscObject)v,&viscuda,VECSEQCUDA,VECMPICUDA,"")); 136*f17f7ee2SSatish Balay PetscCall(PetscObjectTypeCompareAny((PetscObject)s,&siscuda,VECSEQCUDA,VECMPICUDA,"")); 137*f17f7ee2SSatish Balay viscuda = (PetscBool)(viscuda && !v->boundtocpu); 138*f17f7ee2SSatish Balay siscuda = (PetscBool)(siscuda && !s->boundtocpu); 139*f17f7ee2SSatish Balay if (viscuda && siscuda) { 140*f17f7ee2SSatish Balay PetscCall(VecCUDAGetArrayRead(v,&av)); 141*f17f7ee2SSatish Balay PetscCall(VecCUDAGetArrayWrite(s,&as)); 142*f17f7ee2SSatish Balay SignVector_Functor sign_vector(av, as); 143*f17f7ee2SSatish Balay thrust::for_each(thrust::device,thrust::counting_iterator<PetscInt>(0), 144*f17f7ee2SSatish Balay thrust::counting_iterator<PetscInt>(n), sign_vector); 145*f17f7ee2SSatish Balay PetscCall(VecCUDARestoreArrayWrite(s,&as)); 146*f17f7ee2SSatish Balay PetscCall(VecCUDARestoreArrayRead(v,&av)); 147*f17f7ee2SSatish Balay } else 148*f17f7ee2SSatish Balay #endif 149*f17f7ee2SSatish Balay { 150*f17f7ee2SSatish Balay PetscCall(VecGetArrayRead(v,&av)); 151*f17f7ee2SSatish Balay PetscCall(VecGetArrayWrite(s,&as)); 152*f17f7ee2SSatish Balay for (i=0;i<n;i++) as[i] = PetscAbsScalar(av[i]) < 0 ? -1. : 1.; 153*f17f7ee2SSatish Balay PetscCall(VecRestoreArrayWrite(s,&as)); 154*f17f7ee2SSatish Balay PetscCall(VecRestoreArrayRead(v,&av)); 155*f17f7ee2SSatish Balay } 156*f17f7ee2SSatish Balay PetscFunctionReturn(0); 157*f17f7ee2SSatish Balay } 158*f17f7ee2SSatish Balay 159*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 160*f17f7ee2SSatish Balay struct StandardBasis_Functor 161*f17f7ee2SSatish Balay { 162*f17f7ee2SSatish Balay PetscScalar *v; 163*f17f7ee2SSatish Balay PetscInt j; 164*f17f7ee2SSatish Balay 165*f17f7ee2SSatish Balay StandardBasis_Functor(PetscScalar *_v, PetscInt _j) : v(_v), j(_j) {} 166*f17f7ee2SSatish Balay __host__ __device__ void operator()(PetscInt i) 167*f17f7ee2SSatish Balay { 168*f17f7ee2SSatish Balay v[i] = (i == j ? 1 : 0); 169*f17f7ee2SSatish Balay } 170*f17f7ee2SSatish Balay }; 171*f17f7ee2SSatish Balay #endif 172*f17f7ee2SSatish Balay 173*f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSetDelta(Vec x, PetscInt i) 174*f17f7ee2SSatish Balay { 175*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 176*f17f7ee2SSatish Balay PetscBool iscuda; 177*f17f7ee2SSatish Balay #endif 178*f17f7ee2SSatish Balay PetscInt st,en; 179*f17f7ee2SSatish Balay 180*f17f7ee2SSatish Balay PetscFunctionBegin; 181*f17f7ee2SSatish Balay PetscCall(VecGetOwnershipRange(x,&st,&en)); 182*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_CUDA) 183*f17f7ee2SSatish Balay PetscCall(PetscObjectTypeCompareAny((PetscObject)x,&iscuda,VECSEQCUDA,VECMPICUDA,"")); 184*f17f7ee2SSatish Balay iscuda = (PetscBool)(iscuda && !x->boundtocpu); 185*f17f7ee2SSatish Balay if (iscuda) { 186*f17f7ee2SSatish Balay PetscScalar *ax; 187*f17f7ee2SSatish Balay PetscCall(VecCUDAGetArrayWrite(x,&ax)); 188*f17f7ee2SSatish Balay StandardBasis_Functor delta(ax, i-st); 189*f17f7ee2SSatish Balay thrust::for_each(thrust::device,thrust::counting_iterator<PetscInt>(0), 190*f17f7ee2SSatish Balay thrust::counting_iterator<PetscInt>(en-st), delta); 191*f17f7ee2SSatish Balay PetscCall(VecCUDARestoreArrayWrite(x,&ax)); 192*f17f7ee2SSatish Balay } else 193*f17f7ee2SSatish Balay #endif 194*f17f7ee2SSatish Balay { 195*f17f7ee2SSatish Balay PetscCall(VecSet(x,0.)); 196*f17f7ee2SSatish Balay if (st <= i && i < en) { 197*f17f7ee2SSatish Balay PetscCall(VecSetValue(x,i,1.0,INSERT_VALUES)); 198*f17f7ee2SSatish Balay } 199*f17f7ee2SSatish Balay PetscCall(VecAssemblyBegin(x)); 200*f17f7ee2SSatish Balay PetscCall(VecAssemblyEnd(x)); 201*f17f7ee2SSatish Balay } 202*f17f7ee2SSatish Balay PetscFunctionReturn(0); 203*f17f7ee2SSatish Balay } 204*f17f7ee2SSatish Balay 205*f17f7ee2SSatish Balay /* these are approximate norms */ 206*f17f7ee2SSatish Balay /* NORM_2: Estimating the matrix p-norm Nicholas J. Higham 207*f17f7ee2SSatish Balay NORM_1/NORM_INFINITY: A block algorithm for matrix 1-norm estimation, with an application to 1-norm pseudospectra Higham, Nicholas J. and Tisseur, Francoise */ 208*f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode MatApproximateNorm_Private(Mat A, NormType normtype, PetscInt normsamples, PetscReal* n) 209*f17f7ee2SSatish Balay { 210*f17f7ee2SSatish Balay Vec x,y,w,z; 211*f17f7ee2SSatish Balay PetscReal normz,adot; 212*f17f7ee2SSatish Balay PetscScalar dot; 213*f17f7ee2SSatish Balay PetscInt i,j,N,jold = -1; 214*f17f7ee2SSatish Balay PetscBool boundtocpu = PETSC_TRUE; 215*f17f7ee2SSatish Balay 216*f17f7ee2SSatish Balay PetscFunctionBegin; 217*f17f7ee2SSatish Balay #if defined(PETSC_HAVE_DEVICE) 218*f17f7ee2SSatish Balay boundtocpu = A->boundtocpu; 219*f17f7ee2SSatish Balay #endif 220*f17f7ee2SSatish Balay switch (normtype) { 221*f17f7ee2SSatish Balay case NORM_INFINITY: 222*f17f7ee2SSatish Balay case NORM_1: 223*f17f7ee2SSatish Balay if (normsamples < 0) normsamples = 10; /* pure guess */ 224*f17f7ee2SSatish Balay if (normtype == NORM_INFINITY) { 225*f17f7ee2SSatish Balay Mat B; 226*f17f7ee2SSatish Balay PetscCall(MatCreateTranspose(A,&B)); 227*f17f7ee2SSatish Balay A = B; 228*f17f7ee2SSatish Balay } else { 229*f17f7ee2SSatish Balay PetscCall(PetscObjectReference((PetscObject)A)); 230*f17f7ee2SSatish Balay } 231*f17f7ee2SSatish Balay PetscCall(MatCreateVecs(A,&x,&y)); 232*f17f7ee2SSatish Balay PetscCall(MatCreateVecs(A,&z,&w)); 233*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(x,boundtocpu)); 234*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(y,boundtocpu)); 235*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(z,boundtocpu)); 236*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(w,boundtocpu)); 237*f17f7ee2SSatish Balay PetscCall(VecGetSize(x,&N)); 238*f17f7ee2SSatish Balay PetscCall(VecSet(x,1./N)); 239*f17f7ee2SSatish Balay *n = 0.0; 240*f17f7ee2SSatish Balay for (i = 0; i < normsamples; i++) { 241*f17f7ee2SSatish Balay PetscCall(MatMult(A,x,y)); 242*f17f7ee2SSatish Balay PetscCall(VecSign(y,w)); 243*f17f7ee2SSatish Balay PetscCall(MatMultTranspose(A,w,z)); 244*f17f7ee2SSatish Balay PetscCall(VecNorm(z,NORM_INFINITY,&normz)); 245*f17f7ee2SSatish Balay PetscCall(VecDot(x,z,&dot)); 246*f17f7ee2SSatish Balay adot = PetscAbsScalar(dot); 247*f17f7ee2SSatish Balay PetscCall(PetscInfo(A,"%s norm it %" PetscInt_FMT " -> (%g %g)\n",NormTypes[normtype],i,(double)normz,(double)adot)); 248*f17f7ee2SSatish Balay if (normz <= adot && i > 0) { 249*f17f7ee2SSatish Balay PetscCall(VecNorm(y,NORM_1,n)); 250*f17f7ee2SSatish Balay break; 251*f17f7ee2SSatish Balay } 252*f17f7ee2SSatish Balay PetscCall(VecMax(z,&j,&normz)); 253*f17f7ee2SSatish Balay if (j == jold) { 254*f17f7ee2SSatish Balay PetscCall(VecNorm(y,NORM_1,n)); 255*f17f7ee2SSatish Balay PetscCall(PetscInfo(A,"%s norm it %" PetscInt_FMT " -> breakdown (j==jold)\n",NormTypes[normtype],i)); 256*f17f7ee2SSatish Balay break; 257*f17f7ee2SSatish Balay } 258*f17f7ee2SSatish Balay jold = j; 259*f17f7ee2SSatish Balay PetscCall(VecSetDelta(x,j)); 260*f17f7ee2SSatish Balay } 261*f17f7ee2SSatish Balay PetscCall(MatDestroy(&A)); 262*f17f7ee2SSatish Balay PetscCall(VecDestroy(&x)); 263*f17f7ee2SSatish Balay PetscCall(VecDestroy(&w)); 264*f17f7ee2SSatish Balay PetscCall(VecDestroy(&y)); 265*f17f7ee2SSatish Balay PetscCall(VecDestroy(&z)); 266*f17f7ee2SSatish Balay break; 267*f17f7ee2SSatish Balay case NORM_2: 268*f17f7ee2SSatish Balay if (normsamples < 0) normsamples = 20; /* pure guess */ 269*f17f7ee2SSatish Balay PetscCall(MatCreateVecs(A,&x,&y)); 270*f17f7ee2SSatish Balay PetscCall(MatCreateVecs(A,&z,NULL)); 271*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(x,boundtocpu)); 272*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(y,boundtocpu)); 273*f17f7ee2SSatish Balay PetscCall(VecBindToCPU(z,boundtocpu)); 274*f17f7ee2SSatish Balay PetscCall(VecSetRandom(x,NULL)); 275*f17f7ee2SSatish Balay PetscCall(VecNormalize(x,NULL)); 276*f17f7ee2SSatish Balay *n = 0.0; 277*f17f7ee2SSatish Balay for (i = 0; i < normsamples; i++) { 278*f17f7ee2SSatish Balay PetscCall(MatMult(A,x,y)); 279*f17f7ee2SSatish Balay PetscCall(VecNormalize(y,n)); 280*f17f7ee2SSatish Balay PetscCall(MatMultTranspose(A,y,z)); 281*f17f7ee2SSatish Balay PetscCall(VecNorm(z,NORM_2,&normz)); 282*f17f7ee2SSatish Balay PetscCall(VecDot(x,z,&dot)); 283*f17f7ee2SSatish Balay adot = PetscAbsScalar(dot); 284*f17f7ee2SSatish Balay PetscCall(PetscInfo(A,"%s norm it %" PetscInt_FMT " -> %g (%g %g)\n",NormTypes[normtype],i,(double)*n,(double)normz,(double)adot)); 285*f17f7ee2SSatish Balay if (normz <= adot) break; 286*f17f7ee2SSatish Balay if (i < normsamples - 1) { 287*f17f7ee2SSatish Balay Vec t; 288*f17f7ee2SSatish Balay 289*f17f7ee2SSatish Balay PetscCall(VecNormalize(z,NULL)); 290*f17f7ee2SSatish Balay t = x; 291*f17f7ee2SSatish Balay x = z; 292*f17f7ee2SSatish Balay z = t; 293*f17f7ee2SSatish Balay } 294*f17f7ee2SSatish Balay } 295*f17f7ee2SSatish Balay PetscCall(VecDestroy(&x)); 296*f17f7ee2SSatish Balay PetscCall(VecDestroy(&y)); 297*f17f7ee2SSatish Balay PetscCall(VecDestroy(&z)); 298*f17f7ee2SSatish Balay break; 299*f17f7ee2SSatish Balay default: 300*f17f7ee2SSatish Balay SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"%s norm not supported",NormTypes[normtype]); 301*f17f7ee2SSatish Balay } 302*f17f7ee2SSatish Balay PetscCall(PetscInfo(A,"%s norm %g computed in %" PetscInt_FMT " iterations\n",NormTypes[normtype],(double)*n,i)); 303*f17f7ee2SSatish Balay PetscFunctionReturn(0); 304*f17f7ee2SSatish Balay } 305